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)
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 rank %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 /* Currently all group scheme kernels always calculate (shift-)forces */
266 if (flags & GMX_FORCE_FORCES)
268 donb_flags |= GMX_NONBONDED_DO_FORCE;
270 if (flags & GMX_FORCE_VIRIAL)
272 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
274 if (flags & GMX_FORCE_ENERGY)
276 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
278 if (flags & GMX_FORCE_DO_LR)
280 donb_flags |= GMX_NONBONDED_DO_LR;
283 wallcycle_sub_start(wcycle, ewcsNONBONDED);
284 do_nonbonded(fr, x, f, f_longrange, md, excl,
286 lambda, dvdl_nb, -1, -1, donb_flags);
288 /* If we do foreign lambda and we have soft-core interactions
289 * we have to recalculate the (non-linear) energies contributions.
291 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
293 for (i = 0; i < enerd->n_lambda; i++)
295 for (j = 0; j < efptNR; j++)
297 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
299 reset_foreign_enerdata(enerd);
300 do_nonbonded(fr, x, f, f_longrange, md, excl,
301 &(enerd->foreign_grpp), nrnb,
302 lam_i, dvdl_dum, -1, -1,
303 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
304 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
305 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
308 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
312 /* If we are doing GB, calculate bonded forces and apply corrections
313 * to the solvation forces */
314 /* MRS: Eventually, many need to include free energy contribution here! */
315 if (ir->implicit_solvent)
317 wallcycle_sub_start(wcycle, ewcsBONDED);
318 calc_gb_forces(cr, md, born, top, x, f, fr, idef,
319 ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
320 wallcycle_sub_stop(wcycle, ewcsBONDED);
331 if (fepvals->sc_alpha != 0)
333 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
337 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
340 if (fepvals->sc_alpha != 0)
342 /* even though coulomb part is linear, we already added it, beacuse we
343 need to go through the vdw calculation anyway */
345 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
349 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
354 real V_short_range = 0;
355 real dvdl_short_range = 0;
357 for (i = 0; i < enerd->grpp.nener; i++)
361 enerd->grpp.ener[egBHAMSR][i] :
362 enerd->grpp.ener[egLJSR][i])
363 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
365 dvdl_short_range = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
366 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",
375 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
378 /* Shift the coordinates. Must be done before bonded forces and PPPM,
379 * but is also necessary for SHAKE and update, therefore it can NOT
380 * go when no bonded forces have to be evaluated.
383 /* Here sometimes we would not need to shift with NBFonly,
384 * but we do so anyhow for consistency of the returned coordinates.
388 shift_self(graph, box, x);
391 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
395 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
398 /* Check whether we need to do bondeds or correct for exclusions */
400 ((flags & GMX_FORCE_BONDED)
401 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype)))
403 /* Since all atoms are in the rectangular or triclinic unit-cell,
404 * only single box vector shifts (2 in x) are required.
406 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
410 if (flags & GMX_FORCE_BONDED)
412 wallcycle_sub_start(wcycle, ewcsBONDED);
413 calc_bonds(fplog, cr->ms,
414 idef, x, hist, f, fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
415 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
417 fr->bSepDVDL && do_per_step(step, ir->nstlog), step);
419 /* Check if we have to determine energy differences
420 * at foreign lambda's.
422 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
423 idef->ilsort != ilsortNO_FE)
425 if (idef->ilsort != ilsortFE_SORTED)
427 gmx_incons("The bonded interactions are not sorted for free energy");
429 for (i = 0; i < enerd->n_lambda; i++)
431 reset_foreign_enerdata(enerd);
432 for (j = 0; j < efptNR; j++)
434 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
436 calc_bonds_lambda(fplog, idef, x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
437 fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
438 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
439 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
444 wallcycle_sub_stop(wcycle, ewcsBONDED);
450 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
452 real Vlr = 0, Vcorr = 0;
453 real dvdl_long_range = 0;
456 bSB = (ir->nwall == 2);
460 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
461 box_size[ZZ] *= ir->wall_ewald_zfac;
465 /* Do long-range electrostatics and/or LJ-PME, including related short-range
469 clear_mat(fr->vir_el_recip);
470 clear_mat(fr->vir_lj_recip);
472 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
474 real Vlr_q = 0, Vlr_lj = 0, Vcorr_q = 0, Vcorr_lj = 0;
475 real dvdl_long_range_q = 0, dvdl_long_range_lj = 0;
478 if (EEL_PME_EWALD(fr->eeltype) || EVDW_PME(fr->vdwtype))
480 real dvdl_long_range_correction_q = 0;
481 real dvdl_long_range_correction_lj = 0;
482 /* With the Verlet scheme exclusion forces are calculated
483 * in the non-bonded kernel.
485 /* The TPI molecule does not have exclusions with the rest
486 * of the system and no intra-molecular PME grid
487 * contributions will be calculated in
488 * gmx_pme_calc_energy.
490 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
491 ir->ewald_geometry != eewg3D ||
492 ir->epsilon_surface != 0)
496 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
500 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
503 nthreads = gmx_omp_nthreads_get(emntBonded);
504 #pragma omp parallel for num_threads(nthreads) schedule(static)
505 for (t = 0; t < nthreads; t++)
509 tensor *vir_q, *vir_lj;
510 real *Vcorrt_q, *Vcorrt_lj, *dvdlt_q, *dvdlt_lj;
513 fnv = fr->f_novirsum;
514 vir_q = &fr->vir_el_recip;
515 vir_lj = &fr->vir_lj_recip;
517 Vcorrt_lj = &Vcorr_lj;
518 dvdlt_q = &dvdl_long_range_correction_q;
519 dvdlt_lj = &dvdl_long_range_correction_lj;
524 vir_q = &fr->f_t[t].vir_q;
525 vir_lj = &fr->f_t[t].vir_lj;
526 Vcorrt_q = &fr->f_t[t].Vcorr_q;
527 Vcorrt_lj = &fr->f_t[t].Vcorr_lj;
528 dvdlt_q = &fr->f_t[t].dvdl[efptCOUL];
529 dvdlt_lj = &fr->f_t[t].dvdl[efptVDW];
530 for (i = 0; i < fr->natoms_force; i++)
540 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
543 md->nChargePerturbed ? md->chargeB : NULL,
545 md->nTypePerturbed ? md->sqrt_c6B : NULL,
547 md->nTypePerturbed ? md->sigmaB : NULL,
549 md->nTypePerturbed ? md->sigma3B : NULL,
550 ir->cutoff_scheme != ecutsVERLET,
551 excl, x, bSB ? boxs : box, mu_tot,
554 fnv, *vir_q, *vir_lj,
556 lambda[efptCOUL], lambda[efptVDW],
561 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
562 fr->vir_el_recip, fr->vir_lj_recip,
564 &dvdl_long_range_correction_q,
565 &dvdl_long_range_correction_lj,
568 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
571 if (EEL_PME_EWALD(fr->eeltype) && fr->n_tpi == 0)
573 Vcorr_q += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
574 &dvdl_long_range_correction_q,
578 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.", Vcorr_q, dvdl_long_range_correction_q);
579 PRINT_SEPDVDL("Ewald excl. corr. LJ", Vcorr_lj, dvdl_long_range_correction_lj);
580 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction_q;
581 enerd->dvdl_lin[efptVDW] += dvdl_long_range_correction_lj;
584 if ((EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype)))
586 if (cr->duty & DUTY_PME)
588 /* Do reciprocal PME for Coulomb and/or LJ. */
589 assert(fr->n_tpi >= 0);
590 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
592 pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
593 if (EEL_PME(fr->eeltype))
595 pme_flags |= GMX_PME_DO_COULOMB;
597 if (EVDW_PME(fr->vdwtype))
599 pme_flags |= GMX_PME_DO_LJ;
601 if (flags & GMX_FORCE_FORCES)
603 pme_flags |= GMX_PME_CALC_F;
605 if (flags & GMX_FORCE_VIRIAL)
607 pme_flags |= GMX_PME_CALC_ENER_VIR;
611 /* We don't calculate f, but we do want the potential */
612 pme_flags |= GMX_PME_CALC_POT;
614 wallcycle_start(wcycle, ewcPMEMESH);
615 status = gmx_pme_do(fr->pmedata,
616 0, md->homenr - fr->n_tpi,
618 md->chargeA, md->chargeB,
619 md->sqrt_c6A, md->sqrt_c6B,
620 md->sigmaA, md->sigmaB,
621 bSB ? boxs : box, cr,
622 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
623 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
625 fr->vir_el_recip, fr->ewaldcoeff_q,
626 fr->vir_lj_recip, fr->ewaldcoeff_lj,
628 lambda[efptCOUL], lambda[efptVDW],
629 &dvdl_long_range_q, &dvdl_long_range_lj, pme_flags);
630 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
633 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
635 /* We should try to do as little computation after
636 * this as possible, because parallel PME synchronizes
637 * the nodes, so we want all load imbalance of the
638 * rest of the force calculation to be before the PME
639 * call. DD load balancing is done on the whole time
640 * of the force call (without PME).
645 if (EVDW_PME(ir->vdwtype))
648 gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
650 /* Determine the PME grid energy of the test molecule
651 * with the PME grid potential of the other charges.
653 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
654 x + md->homenr - fr->n_tpi,
655 md->chargeA + md->homenr - fr->n_tpi,
658 PRINT_SEPDVDL("PME mesh", Vlr_q + Vlr_lj, dvdl_long_range_q+dvdl_long_range_lj);
662 if (!EEL_PME(fr->eeltype) && EEL_PME_EWALD(fr->eeltype))
664 Vlr_q = do_ewald(ir, x, fr->f_novirsum,
665 md->chargeA, md->chargeB,
666 box_size, cr, md->homenr,
667 fr->vir_el_recip, fr->ewaldcoeff_q,
668 lambda[efptCOUL], &dvdl_long_range_q, fr->ewald_table);
669 PRINT_SEPDVDL("Ewald long-range", Vlr_q, dvdl_long_range_q);
672 /* Note that with separate PME nodes we get the real energies later */
673 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_q;
674 enerd->dvdl_lin[efptVDW] += dvdl_long_range_lj;
675 enerd->term[F_COUL_RECIP] = Vlr_q + Vcorr_q;
676 enerd->term[F_LJ_RECIP] = Vlr_lj + Vcorr_lj;
679 fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
680 Vlr_q, Vcorr_q, enerd->term[F_COUL_RECIP]);
681 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
682 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
683 fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
684 Vlr_lj, Vcorr_lj, enerd->term[F_LJ_RECIP]);
685 pr_rvecs(debug, 0, "vir_lj_recip after corr", fr->vir_lj_recip, DIM);
690 /* Is there a reaction-field exclusion correction needed? */
691 if (EEL_RF(fr->eeltype) && eelRF_NEC != fr->eeltype)
693 /* With the Verlet scheme, exclusion forces are calculated
694 * in the non-bonded kernel.
696 if (ir->cutoff_scheme != ecutsVERLET)
698 real dvdl_rf_excl = 0;
699 enerd->term[F_RF_EXCL] =
700 RF_excl_correction(fr, graph, md, excl, x, f,
701 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
703 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
704 PRINT_SEPDVDL("RF exclusion correction",
705 enerd->term[F_RF_EXCL], dvdl_rf_excl);
714 print_nrnb(debug, nrnb);
722 MPI_Barrier(cr->mpi_comm_mygroup);
725 if (fr->timesteps == 11)
727 fprintf(stderr, "* PP load balancing info: rank %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
728 cr->nodeid, gmx_step_str(fr->timesteps, buf),
729 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
730 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
738 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
743 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
747 for (i = 0; i < F_NRE; i++)
750 enerd->foreign_term[i] = 0;
754 for (i = 0; i < efptNR; i++)
756 enerd->dvdl_lin[i] = 0;
757 enerd->dvdl_nonlin[i] = 0;
763 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
765 enerd->grpp.nener = n2;
766 enerd->foreign_grpp.nener = n2;
767 for (i = 0; (i < egNR); i++)
769 snew(enerd->grpp.ener[i], n2);
770 snew(enerd->foreign_grpp.ener[i], n2);
775 enerd->n_lambda = 1 + n_lambda;
776 snew(enerd->enerpart_lambda, enerd->n_lambda);
784 void destroy_enerdata(gmx_enerdata_t *enerd)
788 for (i = 0; (i < egNR); i++)
790 sfree(enerd->grpp.ener[i]);
793 for (i = 0; (i < egNR); i++)
795 sfree(enerd->foreign_grpp.ener[i]);
800 sfree(enerd->enerpart_lambda);
804 static real sum_v(int n, real v[])
810 for (i = 0; (i < n); i++)
818 void sum_epot(gmx_grppairener_t *grpp, real *epot)
822 /* Accumulate energies */
823 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
824 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
825 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
826 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
827 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
828 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
829 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
830 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
832 /* lattice part of LR doesnt belong to any group
833 * and has been added earlier
835 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
836 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
839 for (i = 0; (i < F_EPOT); i++)
841 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
843 epot[F_EPOT] += epot[i];
848 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
853 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
854 enerd->term[F_DVDL] = 0.0;
855 for (i = 0; i < efptNR; i++)
857 if (fepvals->separate_dvdl[i])
859 /* could this be done more readably/compactly? */
872 index = F_DVDL_BONDED;
874 case (efptRESTRAINT):
875 index = F_DVDL_RESTRAINT;
881 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
884 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
885 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
890 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
893 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
894 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
899 /* Notes on the foreign lambda free energy difference evaluation:
900 * Adding the potential and ekin terms that depend linearly on lambda
901 * as delta lam * dvdl to the energy differences is exact.
902 * For the constraints this is not exact, but we have no other option
903 * without literally changing the lengths and reevaluating the energies at each step.
904 * (try to remedy this post 4.6 - MRS)
905 * For the non-bonded LR term we assume that the soft-core (if present)
906 * no longer affects the energy beyond the short-range cut-off,
907 * which is a very good approximation (except for exotic settings).
908 * (investigate how to overcome this post 4.6 - MRS)
910 if (fepvals->separate_dvdl[efptBONDED])
912 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
916 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
918 enerd->term[F_DVDL_CONSTR] = 0;
920 for (i = 0; i < fepvals->n_lambda; i++)
922 /* note we are iterating over fepvals here!
923 For the current lam, dlam = 0 automatically,
924 so we don't need to add anything to the
925 enerd->enerpart_lambda[0] */
927 /* we don't need to worry about dvdl_lin contributions to dE at
928 current lambda, because the contributions to the current
929 lambda are automatically zeroed */
931 for (j = 0; j < efptNR; j++)
933 /* Note that this loop is over all dhdl components, not just the separated ones */
934 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
935 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
938 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
939 fepvals->all_lambda[j][i], efpt_names[j],
940 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
941 dlam, enerd->dvdl_lin[j]);
948 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
952 /* First reset all foreign energy components. Foreign energies always called on
953 neighbor search steps */
954 for (i = 0; (i < egNR); i++)
956 for (j = 0; (j < enerd->grpp.nener); j++)
958 enerd->foreign_grpp.ener[i][j] = 0.0;
962 /* potential energy components */
963 for (i = 0; (i <= F_EPOT); i++)
965 enerd->foreign_term[i] = 0.0;
969 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
970 gmx_enerdata_t *enerd,
976 /* First reset all energy components, except for the long range terms
977 * on the master at non neighbor search steps, since the long range
978 * terms have already been summed at the last neighbor search step.
980 bKeepLR = (fr->bTwinRange && !bNS);
981 for (i = 0; (i < egNR); i++)
983 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
985 for (j = 0; (j < enerd->grpp.nener); j++)
987 enerd->grpp.ener[i][j] = 0.0;
991 for (i = 0; i < efptNR; i++)
993 enerd->dvdl_lin[i] = 0.0;
994 enerd->dvdl_nonlin[i] = 0.0;
997 /* Normal potential energy components */
998 for (i = 0; (i <= F_EPOT); i++)
1000 enerd->term[i] = 0.0;
1002 /* Initialize the dVdlambda term with the long range contribution */
1003 /* Initialize the dvdl term with the long range contribution */
1004 enerd->term[F_DVDL] = 0.0;
1005 enerd->term[F_DVDL_COUL] = 0.0;
1006 enerd->term[F_DVDL_VDW] = 0.0;
1007 enerd->term[F_DVDL_BONDED] = 0.0;
1008 enerd->term[F_DVDL_RESTRAINT] = 0.0;
1009 enerd->term[F_DKDL] = 0.0;
1010 if (enerd->n_lambda > 0)
1012 for (i = 0; i < enerd->n_lambda; i++)
1014 enerd->enerpart_lambda[i] = 0.0;
1017 /* reset foreign energy data - separate function since we also call it elsewhere */
1018 reset_foreign_enerdata(enerd);