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.
48 #include "nonbonded.h"
60 #include "gmx_omp_nthreads.h"
62 #include "gromacs/legacyheaders/types/commrec.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/pbcutil/ishift.h"
65 #include "gromacs/pbcutil/mshift.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/smalloc.h"
80 gmx_bool bDoLongRangeNS)
86 if (!fr->ns.nblist_initialized)
88 init_neighbor_list(fp, fr, md->homenr);
96 nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
97 bFillGrid, bDoLongRangeNS);
100 fprintf(debug, "nsearch = %d\n", nsearch);
103 /* Check whether we have to do dynamic load balancing */
104 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
105 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
106 &(top->idef),opts->ngener);
108 if (fr->ns.dump_nl > 0)
110 dump_nblist(fp, cr, fr, fr->ns.dump_nl);
114 static void reduce_thread_forces(int n, rvec *f,
115 tensor vir_q, tensor vir_lj,
116 real *Vcorr_q, real *Vcorr_lj,
117 real *dvdl_q, real *dvdl_lj,
118 int nthreads, f_thread_t *f_t)
122 /* This reduction can run over any number of threads */
123 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
124 for (i = 0; i < n; i++)
126 for (t = 1; t < nthreads; t++)
128 rvec_inc(f[i], f_t[t].f[i]);
131 for (t = 1; t < nthreads; t++)
133 *Vcorr_q += f_t[t].Vcorr_q;
134 *Vcorr_lj += f_t[t].Vcorr_lj;
135 *dvdl_q += f_t[t].dvdl[efptCOUL];
136 *dvdl_lj += f_t[t].dvdl[efptVDW];
137 m_add(vir_q, f_t[t].vir_q, vir_q);
138 m_add(vir_lj, f_t[t].vir_lj, vir_lj);
142 void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda)
144 fprintf(fplog, " %-30s V %12.5e dVdl %12.5e\n", s, v, dvdlambda);
147 void do_force_lowlevel(FILE *fplog, gmx_int64_t step,
148 t_forcerec *fr, t_inputrec *ir,
149 t_idef *idef, t_commrec *cr,
150 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
152 rvec x[], history_t *hist,
155 gmx_enerdata_t *enerd,
172 gmx_bool bDoEpot, bSepDVDL, bSB;
178 double clam_i, vlam_i;
179 real dvdl_dum[efptNR], dvdl_nb[efptNR], lam_i[efptNR];
180 real dvdl_q, dvdl_lj;
183 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
186 #define PRINT_SEPDVDL(s, v, dvdlambda) if (bSepDVDL) { gmx_print_sepdvdl(fplog, s, v, dvdlambda); }
188 set_pbc(&pbc, fr->ePBC, box);
190 /* reset free energy components */
191 for (i = 0; i < efptNR; i++)
198 for (i = 0; (i < DIM); i++)
200 box_size[i] = box[i][i];
203 bSepDVDL = (fr->bSepDVDL && do_per_step(step, ir->nstlog));
206 /* do QMMM first if requested */
209 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
214 fprintf(fplog, "Step %s: non-bonded V and dVdl for node %d:\n",
215 gmx_step_str(step, buf), cr->nodeid);
218 /* Call the short range functions all in one go. */
221 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
222 #define TAKETIME FALSE
225 MPI_Barrier(cr->mpi_comm_mygroup);
232 /* foreign lambda component for walls */
233 real dvdl_walls = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
234 enerd->grpp.ener[egLJSR], nrnb);
235 PRINT_SEPDVDL("Walls", 0.0, dvdl_walls);
236 enerd->dvdl_lin[efptVDW] += dvdl_walls;
239 /* If doing GB, reset dvda and calculate the Born radii */
240 if (ir->implicit_solvent)
242 wallcycle_sub_start(wcycle, ewcsNONBONDED);
244 for (i = 0; i < born->nr; i++)
251 calc_gb_rad(cr, fr, ir, top, x, &(fr->gblist), born, md, nrnb);
254 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
258 /* We only do non-bonded calculation with group scheme here, the verlet
259 * calls are done from do_force_cutsVERLET(). */
260 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
263 /* Add short-range interactions */
264 donb_flags |= GMX_NONBONDED_DO_SR;
266 /* Currently all group scheme kernels always calculate (shift-)forces */
267 if (flags & GMX_FORCE_FORCES)
269 donb_flags |= GMX_NONBONDED_DO_FORCE;
271 if (flags & GMX_FORCE_VIRIAL)
273 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
275 if (flags & GMX_FORCE_ENERGY)
277 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
279 if (flags & GMX_FORCE_DO_LR)
281 donb_flags |= GMX_NONBONDED_DO_LR;
284 wallcycle_sub_start(wcycle, ewcsNONBONDED);
285 do_nonbonded(fr, x, f, f_longrange, md, excl,
287 lambda, dvdl_nb, -1, -1, donb_flags);
289 /* If we do foreign lambda and we have soft-core interactions
290 * we have to recalculate the (non-linear) energies contributions.
292 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
294 for (i = 0; i < enerd->n_lambda; i++)
296 for (j = 0; j < efptNR; j++)
298 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
300 reset_foreign_enerdata(enerd);
301 do_nonbonded(fr, x, f, f_longrange, md, excl,
302 &(enerd->foreign_grpp), nrnb,
303 lam_i, dvdl_dum, -1, -1,
304 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
305 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
306 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
309 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
313 /* If we are doing GB, calculate bonded forces and apply corrections
314 * to the solvation forces */
315 /* MRS: Eventually, many need to include free energy contribution here! */
316 if (ir->implicit_solvent)
318 wallcycle_sub_start(wcycle, ewcsBONDED);
319 calc_gb_forces(cr, md, born, top, x, f, fr, idef,
320 ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
321 wallcycle_sub_stop(wcycle, ewcsBONDED);
332 if (fepvals->sc_alpha != 0)
334 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
338 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
341 if (fepvals->sc_alpha != 0)
343 /* even though coulomb part is linear, we already added it, beacuse we
344 need to go through the vdw calculation anyway */
346 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
350 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
355 real V_short_range = 0;
356 real dvdl_short_range = 0;
358 for (i = 0; i < enerd->grpp.nener; i++)
362 enerd->grpp.ener[egBHAMSR][i] :
363 enerd->grpp.ener[egLJSR][i])
364 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
366 dvdl_short_range = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
367 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",
376 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
379 /* Shift the coordinates. Must be done before bonded forces and PPPM,
380 * but is also necessary for SHAKE and update, therefore it can NOT
381 * go when no bonded forces have to be evaluated.
384 /* Here sometimes we would not need to shift with NBFonly,
385 * but we do so anyhow for consistency of the returned coordinates.
389 shift_self(graph, box, x);
392 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
396 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
399 /* Check whether we need to do bondeds or correct for exclusions */
401 ((flags & GMX_FORCE_BONDED)
402 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype)))
404 /* Since all atoms are in the rectangular or triclinic unit-cell,
405 * only single box vector shifts (2 in x) are required.
407 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
411 if (flags & GMX_FORCE_BONDED)
413 wallcycle_sub_start(wcycle, ewcsBONDED);
414 calc_bonds(fplog, cr->ms,
415 idef, x, hist, f, fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
416 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
418 fr->bSepDVDL && do_per_step(step, ir->nstlog), step);
420 /* Check if we have to determine energy differences
421 * at foreign lambda's.
423 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
424 idef->ilsort != ilsortNO_FE)
426 if (idef->ilsort != ilsortFE_SORTED)
428 gmx_incons("The bonded interactions are not sorted for free energy");
430 for (i = 0; i < enerd->n_lambda; i++)
432 reset_foreign_enerdata(enerd);
433 for (j = 0; j < efptNR; j++)
435 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
437 calc_bonds_lambda(fplog, idef, x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
438 fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
439 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
440 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
445 wallcycle_sub_stop(wcycle, ewcsBONDED);
451 clear_mat(fr->vir_el_recip);
452 clear_mat(fr->vir_lj_recip);
454 /* Do long-range electrostatics and/or LJ-PME, including related short-range
457 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
459 real Vlr = 0, Vcorr = 0;
460 real dvdl_long_range = 0;
462 real Vlr_q = 0, Vlr_lj = 0, Vcorr_q = 0, Vcorr_lj = 0;
463 real dvdl_long_range_q = 0, dvdl_long_range_lj = 0;
465 bSB = (ir->nwall == 2);
469 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
470 box_size[ZZ] *= ir->wall_ewald_zfac;
473 if (EEL_PME_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++)
535 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
538 md->nChargePerturbed ? md->chargeB : NULL,
540 md->nTypePerturbed ? md->sqrt_c6B : NULL,
542 md->nTypePerturbed ? md->sigmaB : NULL,
544 md->nTypePerturbed ? md->sigma3B : NULL,
545 ir->cutoff_scheme != ecutsVERLET,
546 excl, x, bSB ? boxs : box, mu_tot,
549 fnv, *vir_q, *vir_lj,
551 lambda[efptCOUL], lambda[efptVDW],
556 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
557 fr->vir_el_recip, fr->vir_lj_recip,
559 &dvdl_long_range_correction_q,
560 &dvdl_long_range_correction_lj,
563 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
566 if (EEL_PME_EWALD(fr->eeltype) && fr->n_tpi == 0)
568 Vcorr_q += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
569 &dvdl_long_range_correction_q,
573 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr_q, dvdl_long_range_correction_q);
574 PRINT_SEPDVDL("Ewald excl. corr. LJ", Vcorr_lj, dvdl_long_range_correction_lj);
575 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction_q;
576 enerd->dvdl_lin[efptVDW] += dvdl_long_range_correction_lj;
578 if ((EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype)) && (cr->duty & DUTY_PME))
580 /* Do reciprocal PME for Coulomb and/or LJ. */
581 assert(fr->n_tpi >= 0);
582 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
584 pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
585 if (EEL_PME(fr->eeltype))
587 pme_flags |= GMX_PME_DO_COULOMB;
589 if (EVDW_PME(fr->vdwtype))
591 pme_flags |= GMX_PME_DO_LJ;
593 if (flags & GMX_FORCE_FORCES)
595 pme_flags |= GMX_PME_CALC_F;
597 if (flags & GMX_FORCE_VIRIAL)
599 pme_flags |= GMX_PME_CALC_ENER_VIR;
603 /* We don't calculate f, but we do want the potential */
604 pme_flags |= GMX_PME_CALC_POT;
606 wallcycle_start(wcycle, ewcPMEMESH);
607 status = gmx_pme_do(fr->pmedata,
608 0, md->homenr - fr->n_tpi,
610 md->chargeA, md->chargeB,
611 md->sqrt_c6A, md->sqrt_c6B,
612 md->sigmaA, md->sigmaB,
613 bSB ? boxs : box, cr,
614 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
615 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
617 fr->vir_el_recip, fr->ewaldcoeff_q,
618 fr->vir_lj_recip, fr->ewaldcoeff_lj,
620 lambda[efptCOUL], lambda[efptVDW],
621 &dvdl_long_range_q, &dvdl_long_range_lj, pme_flags);
622 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
625 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
627 /* We should try to do as little computation after
628 * this as possible, because parallel PME synchronizes
629 * the nodes, so we want all load imbalance of the
630 * rest of the force calculation to be before the PME
631 * call. DD load balancing is done on the whole time
632 * of the force call (without PME).
637 if (EVDW_PME(ir->vdwtype))
640 gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
642 /* Determine the PME grid energy of the test molecule
643 * with the PME grid potential of the other charges.
645 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
646 x + md->homenr - fr->n_tpi,
647 md->chargeA + md->homenr - fr->n_tpi,
650 PRINT_SEPDVDL("PME mesh", Vlr_q + Vlr_lj, dvdl_long_range_q+dvdl_long_range_lj);
654 if (!EEL_PME(fr->eeltype) && EEL_PME_EWALD(fr->eeltype))
656 Vlr_q = do_ewald(ir, x, fr->f_novirsum,
657 md->chargeA, md->chargeB,
658 box_size, cr, md->homenr,
659 fr->vir_el_recip, fr->ewaldcoeff_q,
660 lambda[efptCOUL], &dvdl_long_range_q, fr->ewald_table);
661 PRINT_SEPDVDL("Ewald long-range", Vlr_q, dvdl_long_range_q);
664 /* Note that with separate PME nodes we get the real energies later */
665 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_q;
666 enerd->dvdl_lin[efptVDW] += dvdl_long_range_lj;
667 enerd->term[F_COUL_RECIP] = Vlr_q + Vcorr_q;
668 enerd->term[F_LJ_RECIP] = Vlr_lj + Vcorr_lj;
671 fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
672 Vlr_q, Vcorr_q, enerd->term[F_COUL_RECIP]);
673 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
674 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
675 fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
676 Vlr_lj, Vcorr_lj, enerd->term[F_LJ_RECIP]);
677 pr_rvecs(debug, 0, "vir_lj_recip after corr", fr->vir_lj_recip, DIM);
682 /* Is there a reaction-field exclusion correction needed? */
683 if (EEL_RF(fr->eeltype) && eelRF_NEC != fr->eeltype)
685 /* With the Verlet scheme, exclusion forces are calculated
686 * in the non-bonded kernel.
688 if (ir->cutoff_scheme != ecutsVERLET)
690 real dvdl_rf_excl = 0;
691 enerd->term[F_RF_EXCL] =
692 RF_excl_correction(fr, graph, md, excl, x, f,
693 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
695 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
696 PRINT_SEPDVDL("RF exclusion correction",
697 enerd->term[F_RF_EXCL], dvdl_rf_excl);
706 print_nrnb(debug, nrnb);
714 MPI_Barrier(cr->mpi_comm_mygroup);
717 if (fr->timesteps == 11)
719 fprintf(stderr, "* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
720 cr->nodeid, gmx_step_str(fr->timesteps, buf),
721 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
722 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
730 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
735 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
739 for (i = 0; i < F_NRE; i++)
742 enerd->foreign_term[i] = 0;
746 for (i = 0; i < efptNR; i++)
748 enerd->dvdl_lin[i] = 0;
749 enerd->dvdl_nonlin[i] = 0;
755 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
757 enerd->grpp.nener = n2;
758 enerd->foreign_grpp.nener = n2;
759 for (i = 0; (i < egNR); i++)
761 snew(enerd->grpp.ener[i], n2);
762 snew(enerd->foreign_grpp.ener[i], n2);
767 enerd->n_lambda = 1 + n_lambda;
768 snew(enerd->enerpart_lambda, enerd->n_lambda);
776 void destroy_enerdata(gmx_enerdata_t *enerd)
780 for (i = 0; (i < egNR); i++)
782 sfree(enerd->grpp.ener[i]);
785 for (i = 0; (i < egNR); i++)
787 sfree(enerd->foreign_grpp.ener[i]);
792 sfree(enerd->enerpart_lambda);
796 static real sum_v(int n, real v[])
802 for (i = 0; (i < n); i++)
810 void sum_epot(gmx_grppairener_t *grpp, real *epot)
814 /* Accumulate energies */
815 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
816 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
817 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
818 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
819 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
820 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
821 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
822 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
824 /* lattice part of LR doesnt belong to any group
825 * and has been added earlier
827 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
828 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
831 for (i = 0; (i < F_EPOT); i++)
833 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
835 epot[F_EPOT] += epot[i];
840 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
845 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
846 enerd->term[F_DVDL] = 0.0;
847 for (i = 0; i < efptNR; i++)
849 if (fepvals->separate_dvdl[i])
851 /* could this be done more readably/compactly? */
864 index = F_DVDL_BONDED;
866 case (efptRESTRAINT):
867 index = F_DVDL_RESTRAINT;
873 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
876 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
877 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
882 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
885 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
886 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
891 /* Notes on the foreign lambda free energy difference evaluation:
892 * Adding the potential and ekin terms that depend linearly on lambda
893 * as delta lam * dvdl to the energy differences is exact.
894 * For the constraints this is not exact, but we have no other option
895 * without literally changing the lengths and reevaluating the energies at each step.
896 * (try to remedy this post 4.6 - MRS)
897 * For the non-bonded LR term we assume that the soft-core (if present)
898 * no longer affects the energy beyond the short-range cut-off,
899 * which is a very good approximation (except for exotic settings).
900 * (investigate how to overcome this post 4.6 - MRS)
902 if (fepvals->separate_dvdl[efptBONDED])
904 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
908 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
910 enerd->term[F_DVDL_CONSTR] = 0;
912 for (i = 0; i < fepvals->n_lambda; i++)
914 /* note we are iterating over fepvals here!
915 For the current lam, dlam = 0 automatically,
916 so we don't need to add anything to the
917 enerd->enerpart_lambda[0] */
919 /* we don't need to worry about dvdl_lin contributions to dE at
920 current lambda, because the contributions to the current
921 lambda are automatically zeroed */
923 for (j = 0; j < efptNR; j++)
925 /* Note that this loop is over all dhdl components, not just the separated ones */
926 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
927 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
930 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
931 fepvals->all_lambda[j][i], efpt_names[j],
932 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
933 dlam, enerd->dvdl_lin[j]);
940 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
944 /* First reset all foreign energy components. Foreign energies always called on
945 neighbor search steps */
946 for (i = 0; (i < egNR); i++)
948 for (j = 0; (j < enerd->grpp.nener); j++)
950 enerd->foreign_grpp.ener[i][j] = 0.0;
954 /* potential energy components */
955 for (i = 0; (i <= F_EPOT); i++)
957 enerd->foreign_term[i] = 0.0;
961 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
962 gmx_enerdata_t *enerd,
968 /* First reset all energy components, except for the long range terms
969 * on the master at non neighbor search steps, since the long range
970 * terms have already been summed at the last neighbor search step.
972 bKeepLR = (fr->bTwinRange && !bNS);
973 for (i = 0; (i < egNR); i++)
975 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
977 for (j = 0; (j < enerd->grpp.nener); j++)
979 enerd->grpp.ener[i][j] = 0.0;
983 for (i = 0; i < efptNR; i++)
985 enerd->dvdl_lin[i] = 0.0;
986 enerd->dvdl_nonlin[i] = 0.0;
989 /* Normal potential energy components */
990 for (i = 0; (i <= F_EPOT); i++)
992 enerd->term[i] = 0.0;
994 /* Initialize the dVdlambda term with the long range contribution */
995 /* Initialize the dvdl term with the long range contribution */
996 enerd->term[F_DVDL] = 0.0;
997 enerd->term[F_DVDL_COUL] = 0.0;
998 enerd->term[F_DVDL_VDW] = 0.0;
999 enerd->term[F_DVDL_BONDED] = 0.0;
1000 enerd->term[F_DVDL_RESTRAINT] = 0.0;
1001 enerd->term[F_DKDL] = 0.0;
1002 if (enerd->n_lambda > 0)
1004 for (i = 0; i < enerd->n_lambda; i++)
1006 enerd->enerpart_lambda[i] = 0.0;
1009 /* reset foreign energy data - separate function since we also call it elsewhere */
1010 reset_foreign_enerdata(enerd);