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.
47 #include "gromacs/utility/smalloc.h"
51 #include "nonbonded.h"
65 #include "gmx_omp_nthreads.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gmx_fatal.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)
120 int nthreads_loop gmx_unused;
122 /* This reduction can run over any number of threads */
123 nthreads_loop = gmx_omp_nthreads_get(emntBonded);
124 #pragma omp parallel for num_threads(nthreads_loop) private(t) schedule(static)
125 for (i = 0; i < n; i++)
127 for (t = 1; t < nthreads; t++)
129 rvec_inc(f[i], f_t[t].f[i]);
132 for (t = 1; t < nthreads; t++)
134 *Vcorr_q += f_t[t].Vcorr_q;
135 *Vcorr_lj += f_t[t].Vcorr_lj;
136 *dvdl_q += f_t[t].dvdl[efptCOUL];
137 *dvdl_lj += f_t[t].dvdl[efptVDW];
138 m_add(vir_q, f_t[t].vir_q, vir_q);
139 m_add(vir_lj, f_t[t].vir_lj, vir_lj);
143 void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda)
145 fprintf(fplog, " %-30s V %12.5e dVdl %12.5e\n", s, v, dvdlambda);
148 void do_force_lowlevel(FILE *fplog, gmx_int64_t step,
149 t_forcerec *fr, t_inputrec *ir,
150 t_idef *idef, t_commrec *cr,
151 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
153 rvec x[], history_t *hist,
156 gmx_enerdata_t *enerd,
173 gmx_bool bDoEpot, bSepDVDL, bSB;
179 double clam_i, vlam_i;
180 real dvdl_dum[efptNR], dvdl_nb[efptNR], lam_i[efptNR];
181 real dvdl_q, dvdl_lj;
184 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
187 #define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) { gmx_print_sepdvdl(fplog, s, v, dvdlambda); }
189 set_pbc(&pbc, fr->ePBC, box);
191 /* reset free energy components */
192 for (i = 0; i < efptNR; i++)
199 for (i = 0; (i < DIM); i++)
201 box_size[i] = box[i][i];
204 bSepDVDL = (fr->bSepDVDL && do_per_step(step, ir->nstlog));
207 /* do QMMM first if requested */
210 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
215 fprintf(fplog, "Step %s: non-bonded V and dVdl for rank %d:\n",
216 gmx_step_str(step, buf), cr->nodeid);
219 /* Call the short range functions all in one go. */
222 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
223 #define TAKETIME FALSE
226 MPI_Barrier(cr->mpi_comm_mygroup);
233 /* foreign lambda component for walls */
234 real dvdl_walls = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
235 enerd->grpp.ener[egLJSR], nrnb);
236 PRINT_SEPDVDL("Walls", 0.0, dvdl_walls);
237 enerd->dvdl_lin[efptVDW] += dvdl_walls;
240 /* If doing GB, reset dvda and calculate the Born radii */
241 if (ir->implicit_solvent)
243 wallcycle_sub_start(wcycle, ewcsNONBONDED);
245 for (i = 0; i < born->nr; i++)
252 calc_gb_rad(cr, fr, ir, top, x, &(fr->gblist), born, md, nrnb);
255 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
259 /* We only do non-bonded calculation with group scheme here, the verlet
260 * calls are done from do_force_cutsVERLET(). */
261 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
264 /* Add short-range interactions */
265 donb_flags |= GMX_NONBONDED_DO_SR;
267 /* Currently all group scheme kernels always calculate (shift-)forces */
268 if (flags & GMX_FORCE_FORCES)
270 donb_flags |= GMX_NONBONDED_DO_FORCE;
272 if (flags & GMX_FORCE_VIRIAL)
274 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
276 if (flags & GMX_FORCE_ENERGY)
278 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
280 if (flags & GMX_FORCE_DO_LR)
282 donb_flags |= GMX_NONBONDED_DO_LR;
285 wallcycle_sub_start(wcycle, ewcsNONBONDED);
286 do_nonbonded(fr, x, f, f_longrange, md, excl,
288 lambda, dvdl_nb, -1, -1, donb_flags);
290 /* If we do foreign lambda and we have soft-core interactions
291 * we have to recalculate the (non-linear) energies contributions.
293 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
295 for (i = 0; i < enerd->n_lambda; i++)
297 for (j = 0; j < efptNR; j++)
299 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
301 reset_foreign_enerdata(enerd);
302 do_nonbonded(fr, x, f, f_longrange, md, excl,
303 &(enerd->foreign_grpp), nrnb,
304 lam_i, dvdl_dum, -1, -1,
305 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
306 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
307 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
310 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
314 /* If we are doing GB, calculate bonded forces and apply corrections
315 * to the solvation forces */
316 /* MRS: Eventually, many need to include free energy contribution here! */
317 if (ir->implicit_solvent)
319 wallcycle_sub_start(wcycle, ewcsBONDED);
320 calc_gb_forces(cr, md, born, top, x, f, fr, idef,
321 ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
322 wallcycle_sub_stop(wcycle, ewcsBONDED);
333 if (fepvals->sc_alpha != 0)
335 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
339 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
342 if (fepvals->sc_alpha != 0)
344 /* even though coulomb part is linear, we already added it, beacuse we
345 need to go through the vdw calculation anyway */
347 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
351 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
356 real V_short_range = 0;
357 real dvdl_short_range = 0;
359 for (i = 0; i < enerd->grpp.nener; i++)
363 enerd->grpp.ener[egBHAMSR][i] :
364 enerd->grpp.ener[egLJSR][i])
365 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
367 dvdl_short_range = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
368 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",
377 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
380 /* Shift the coordinates. Must be done before bonded forces and PPPM,
381 * but is also necessary for SHAKE and update, therefore it can NOT
382 * go when no bonded forces have to be evaluated.
385 /* Here sometimes we would not need to shift with NBFonly,
386 * but we do so anyhow for consistency of the returned coordinates.
390 shift_self(graph, box, x);
393 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
397 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
400 /* Check whether we need to do bondeds or correct for exclusions */
402 ((flags & GMX_FORCE_BONDED)
403 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype)))
405 /* Since all atoms are in the rectangular or triclinic unit-cell,
406 * only single box vector shifts (2 in x) are required.
408 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
412 if (flags & GMX_FORCE_BONDED)
414 wallcycle_sub_start(wcycle, ewcsBONDED);
415 calc_bonds(fplog, cr->ms,
416 idef, x, hist, f, fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
417 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
419 fr->bSepDVDL && do_per_step(step, ir->nstlog), step);
421 /* Check if we have to determine energy differences
422 * at foreign lambda's.
424 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
425 idef->ilsort != ilsortNO_FE)
427 if (idef->ilsort != ilsortFE_SORTED)
429 gmx_incons("The bonded interactions are not sorted for free energy");
431 for (i = 0; i < enerd->n_lambda; i++)
433 reset_foreign_enerdata(enerd);
434 for (j = 0; j < efptNR; j++)
436 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
438 calc_bonds_lambda(fplog, idef, x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
439 fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
440 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
441 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
446 wallcycle_sub_stop(wcycle, ewcsBONDED);
452 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
454 real Vlr = 0, Vcorr = 0;
455 real dvdl_long_range = 0;
458 bSB = (ir->nwall == 2);
462 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
463 box_size[ZZ] *= ir->wall_ewald_zfac;
467 /* Do long-range electrostatics and/or LJ-PME, including related short-range
471 clear_mat(fr->vir_el_recip);
472 clear_mat(fr->vir_lj_recip);
474 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
476 real Vlr_q = 0, Vlr_lj = 0, Vcorr_q = 0, Vcorr_lj = 0;
477 real dvdl_long_range_q = 0, dvdl_long_range_lj = 0;
480 if (EEL_PME_EWALD(fr->eeltype) || EVDW_PME(fr->vdwtype))
482 real dvdl_long_range_correction_q = 0;
483 real dvdl_long_range_correction_lj = 0;
484 /* With the Verlet scheme exclusion forces are calculated
485 * in the non-bonded kernel.
487 /* The TPI molecule does not have exclusions with the rest
488 * of the system and no intra-molecular PME grid
489 * contributions will be calculated in
490 * gmx_pme_calc_energy.
492 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
493 ir->ewald_geometry != eewg3D ||
494 ir->epsilon_surface != 0)
498 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
502 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
505 nthreads = gmx_omp_nthreads_get(emntBonded);
506 #pragma omp parallel for num_threads(nthreads) schedule(static)
507 for (t = 0; t < nthreads; t++)
511 tensor *vir_q, *vir_lj;
512 real *Vcorrt_q, *Vcorrt_lj, *dvdlt_q, *dvdlt_lj;
515 fnv = fr->f_novirsum;
516 vir_q = &fr->vir_el_recip;
517 vir_lj = &fr->vir_lj_recip;
519 Vcorrt_lj = &Vcorr_lj;
520 dvdlt_q = &dvdl_long_range_correction_q;
521 dvdlt_lj = &dvdl_long_range_correction_lj;
526 vir_q = &fr->f_t[t].vir_q;
527 vir_lj = &fr->f_t[t].vir_lj;
528 Vcorrt_q = &fr->f_t[t].Vcorr_q;
529 Vcorrt_lj = &fr->f_t[t].Vcorr_lj;
530 dvdlt_q = &fr->f_t[t].dvdl[efptCOUL];
531 dvdlt_lj = &fr->f_t[t].dvdl[efptVDW];
532 for (i = 0; i < fr->natoms_force; i++)
542 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
545 md->nChargePerturbed ? md->chargeB : NULL,
547 md->nTypePerturbed ? md->sqrt_c6B : NULL,
549 md->nTypePerturbed ? md->sigmaB : NULL,
551 md->nTypePerturbed ? md->sigma3B : NULL,
552 ir->cutoff_scheme != ecutsVERLET,
553 excl, x, bSB ? boxs : box, mu_tot,
556 fnv, *vir_q, *vir_lj,
558 lambda[efptCOUL], lambda[efptVDW],
563 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
564 fr->vir_el_recip, fr->vir_lj_recip,
566 &dvdl_long_range_correction_q,
567 &dvdl_long_range_correction_lj,
570 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
573 if (EEL_PME_EWALD(fr->eeltype) && fr->n_tpi == 0)
575 Vcorr_q += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
576 &dvdl_long_range_correction_q,
580 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr_q, dvdl_long_range_correction_q);
581 PRINT_SEPDVDL("Ewald excl. corr. LJ", Vcorr_lj, dvdl_long_range_correction_lj);
582 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction_q;
583 enerd->dvdl_lin[efptVDW] += dvdl_long_range_correction_lj;
586 if ((EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype)))
588 if (cr->duty & DUTY_PME)
590 /* Do reciprocal PME for Coulomb and/or LJ. */
591 assert(fr->n_tpi >= 0);
592 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
594 pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
595 if (EEL_PME(fr->eeltype))
597 pme_flags |= GMX_PME_DO_COULOMB;
599 if (EVDW_PME(fr->vdwtype))
601 pme_flags |= GMX_PME_DO_LJ;
603 if (flags & GMX_FORCE_FORCES)
605 pme_flags |= GMX_PME_CALC_F;
607 if (flags & GMX_FORCE_VIRIAL)
609 pme_flags |= GMX_PME_CALC_ENER_VIR;
613 /* We don't calculate f, but we do want the potential */
614 pme_flags |= GMX_PME_CALC_POT;
616 wallcycle_start(wcycle, ewcPMEMESH);
617 status = gmx_pme_do(fr->pmedata,
618 0, md->homenr - fr->n_tpi,
620 md->chargeA, md->chargeB,
621 md->sqrt_c6A, md->sqrt_c6B,
622 md->sigmaA, md->sigmaB,
623 bSB ? boxs : box, cr,
624 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
625 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
627 fr->vir_el_recip, fr->ewaldcoeff_q,
628 fr->vir_lj_recip, fr->ewaldcoeff_lj,
630 lambda[efptCOUL], lambda[efptVDW],
631 &dvdl_long_range_q, &dvdl_long_range_lj, pme_flags);
632 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
635 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
637 /* We should try to do as little computation after
638 * this as possible, because parallel PME synchronizes
639 * the nodes, so we want all load imbalance of the
640 * rest of the force calculation to be before the PME
641 * call. DD load balancing is done on the whole time
642 * of the force call (without PME).
647 if (EVDW_PME(ir->vdwtype))
650 gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
652 /* Determine the PME grid energy of the test molecule
653 * with the PME grid potential of the other charges.
655 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
656 x + md->homenr - fr->n_tpi,
657 md->chargeA + md->homenr - fr->n_tpi,
660 PRINT_SEPDVDL("PME mesh", Vlr_q + Vlr_lj, dvdl_long_range_q+dvdl_long_range_lj);
664 if (!EEL_PME(fr->eeltype) && EEL_PME_EWALD(fr->eeltype))
666 Vlr_q = do_ewald(ir, x, fr->f_novirsum,
667 md->chargeA, md->chargeB,
668 box_size, cr, md->homenr,
669 fr->vir_el_recip, fr->ewaldcoeff_q,
670 lambda[efptCOUL], &dvdl_long_range_q, fr->ewald_table);
671 PRINT_SEPDVDL("Ewald long-range", Vlr_q, dvdl_long_range_q);
674 /* Note that with separate PME nodes we get the real energies later */
675 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_q;
676 enerd->dvdl_lin[efptVDW] += dvdl_long_range_lj;
677 enerd->term[F_COUL_RECIP] = Vlr_q + Vcorr_q;
678 enerd->term[F_LJ_RECIP] = Vlr_lj + Vcorr_lj;
681 fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
682 Vlr_q, Vcorr_q, enerd->term[F_COUL_RECIP]);
683 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
684 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
685 fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
686 Vlr_lj, Vcorr_lj, enerd->term[F_LJ_RECIP]);
687 pr_rvecs(debug, 0, "vir_lj_recip after corr", fr->vir_lj_recip, DIM);
692 /* Is there a reaction-field exclusion correction needed? */
693 if (EEL_RF(fr->eeltype) && eelRF_NEC != fr->eeltype)
695 /* With the Verlet scheme, exclusion forces are calculated
696 * in the non-bonded kernel.
698 if (ir->cutoff_scheme != ecutsVERLET)
700 real dvdl_rf_excl = 0;
701 enerd->term[F_RF_EXCL] =
702 RF_excl_correction(fr, graph, md, excl, x, f,
703 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
705 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
706 PRINT_SEPDVDL("RF exclusion correction",
707 enerd->term[F_RF_EXCL], dvdl_rf_excl);
716 print_nrnb(debug, nrnb);
724 MPI_Barrier(cr->mpi_comm_mygroup);
727 if (fr->timesteps == 11)
729 fprintf(stderr, "* PP load balancing info: rank %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
730 cr->nodeid, gmx_step_str(fr->timesteps, buf),
731 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
732 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
740 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
745 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
749 for (i = 0; i < F_NRE; i++)
752 enerd->foreign_term[i] = 0;
756 for (i = 0; i < efptNR; i++)
758 enerd->dvdl_lin[i] = 0;
759 enerd->dvdl_nonlin[i] = 0;
765 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
767 enerd->grpp.nener = n2;
768 enerd->foreign_grpp.nener = n2;
769 for (i = 0; (i < egNR); i++)
771 snew(enerd->grpp.ener[i], n2);
772 snew(enerd->foreign_grpp.ener[i], n2);
777 enerd->n_lambda = 1 + n_lambda;
778 snew(enerd->enerpart_lambda, enerd->n_lambda);
786 void destroy_enerdata(gmx_enerdata_t *enerd)
790 for (i = 0; (i < egNR); i++)
792 sfree(enerd->grpp.ener[i]);
795 for (i = 0; (i < egNR); i++)
797 sfree(enerd->foreign_grpp.ener[i]);
802 sfree(enerd->enerpart_lambda);
806 static real sum_v(int n, real v[])
812 for (i = 0; (i < n); i++)
820 void sum_epot(gmx_grppairener_t *grpp, real *epot)
824 /* Accumulate energies */
825 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
826 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
827 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
828 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
829 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
830 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
831 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
832 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
834 /* lattice part of LR doesnt belong to any group
835 * and has been added earlier
837 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
838 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
841 for (i = 0; (i < F_EPOT); i++)
843 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
845 epot[F_EPOT] += epot[i];
850 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
855 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
856 enerd->term[F_DVDL] = 0.0;
857 for (i = 0; i < efptNR; i++)
859 if (fepvals->separate_dvdl[i])
861 /* could this be done more readably/compactly? */
874 index = F_DVDL_BONDED;
876 case (efptRESTRAINT):
877 index = F_DVDL_RESTRAINT;
883 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
886 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
887 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
892 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
895 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
896 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
901 /* Notes on the foreign lambda free energy difference evaluation:
902 * Adding the potential and ekin terms that depend linearly on lambda
903 * as delta lam * dvdl to the energy differences is exact.
904 * For the constraints this is not exact, but we have no other option
905 * without literally changing the lengths and reevaluating the energies at each step.
906 * (try to remedy this post 4.6 - MRS)
907 * For the non-bonded LR term we assume that the soft-core (if present)
908 * no longer affects the energy beyond the short-range cut-off,
909 * which is a very good approximation (except for exotic settings).
910 * (investigate how to overcome this post 4.6 - MRS)
912 if (fepvals->separate_dvdl[efptBONDED])
914 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
918 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
920 enerd->term[F_DVDL_CONSTR] = 0;
922 for (i = 0; i < fepvals->n_lambda; i++)
924 /* note we are iterating over fepvals here!
925 For the current lam, dlam = 0 automatically,
926 so we don't need to add anything to the
927 enerd->enerpart_lambda[0] */
929 /* we don't need to worry about dvdl_lin contributions to dE at
930 current lambda, because the contributions to the current
931 lambda are automatically zeroed */
933 for (j = 0; j < efptNR; j++)
935 /* Note that this loop is over all dhdl components, not just the separated ones */
936 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
937 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
940 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
941 fepvals->all_lambda[j][i], efpt_names[j],
942 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
943 dlam, enerd->dvdl_lin[j]);
950 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
954 /* First reset all foreign energy components. Foreign energies always called on
955 neighbor search steps */
956 for (i = 0; (i < egNR); i++)
958 for (j = 0; (j < enerd->grpp.nener); j++)
960 enerd->foreign_grpp.ener[i][j] = 0.0;
964 /* potential energy components */
965 for (i = 0; (i <= F_EPOT); i++)
967 enerd->foreign_term[i] = 0.0;
971 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
972 gmx_enerdata_t *enerd,
978 /* First reset all energy components, except for the long range terms
979 * on the master at non neighbor search steps, since the long range
980 * terms have already been summed at the last neighbor search step.
982 bKeepLR = (fr->bTwinRange && !bNS);
983 for (i = 0; (i < egNR); i++)
985 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
987 for (j = 0; (j < enerd->grpp.nener); j++)
989 enerd->grpp.ener[i][j] = 0.0;
993 for (i = 0; i < efptNR; i++)
995 enerd->dvdl_lin[i] = 0.0;
996 enerd->dvdl_nonlin[i] = 0.0;
999 /* Normal potential energy components */
1000 for (i = 0; (i <= F_EPOT); i++)
1002 enerd->term[i] = 0.0;
1004 /* Initialize the dVdlambda term with the long range contribution */
1005 /* Initialize the dvdl term with the long range contribution */
1006 enerd->term[F_DVDL] = 0.0;
1007 enerd->term[F_DVDL_COUL] = 0.0;
1008 enerd->term[F_DVDL_VDW] = 0.0;
1009 enerd->term[F_DVDL_BONDED] = 0.0;
1010 enerd->term[F_DVDL_RESTRAINT] = 0.0;
1011 enerd->term[F_DKDL] = 0.0;
1012 if (enerd->n_lambda > 0)
1014 for (i = 0; i < enerd->n_lambda; i++)
1016 enerd->enerpart_lambda[i] = 0.0;
1019 /* reset foreign energy data - separate function since we also call it elsewhere */
1020 reset_foreign_enerdata(enerd);