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.
39 #include "gromacs/legacyheaders/force.h"
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/ewald/ewald.h"
49 #include "gromacs/ewald/long-range-correction.h"
50 #include "gromacs/ewald/pme.h"
51 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
52 #include "gromacs/legacyheaders/macros.h"
53 #include "gromacs/legacyheaders/mdrun.h"
54 #include "gromacs/legacyheaders/names.h"
55 #include "gromacs/legacyheaders/network.h"
56 #include "gromacs/legacyheaders/nonbonded.h"
57 #include "gromacs/legacyheaders/nrnb.h"
58 #include "gromacs/legacyheaders/ns.h"
59 #include "gromacs/legacyheaders/qmmm.h"
60 #include "gromacs/legacyheaders/txtdump.h"
61 #include "gromacs/legacyheaders/typedefs.h"
62 #include "gromacs/legacyheaders/types/commrec.h"
63 #include "gromacs/listed-forces/listed-forces.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/pbcutil/ishift.h"
66 #include "gromacs/pbcutil/mshift.h"
67 #include "gromacs/pbcutil/pbc.h"
68 #include "gromacs/timing/wallcycle.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/smalloc.h"
81 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)
121 int nthreads_loop gmx_unused;
123 // cppcheck-suppress unreadVariable
124 nthreads_loop = gmx_omp_nthreads_get(emntBonded);
125 /* This reduction can run over any number of threads */
126 #pragma omp parallel for num_threads(nthreads_loop) private(t) schedule(static)
127 for (i = 0; i < n; i++)
129 for (t = 1; t < nthreads; t++)
131 rvec_inc(f[i], f_t[t].f[i]);
134 for (t = 1; t < nthreads; t++)
136 *Vcorr_q += f_t[t].Vcorr_q;
137 *Vcorr_lj += f_t[t].Vcorr_lj;
138 *dvdl_q += f_t[t].dvdl[efptCOUL];
139 *dvdl_lj += f_t[t].dvdl[efptVDW];
140 m_add(vir_q, f_t[t].vir_q, vir_q);
141 m_add(vir_lj, f_t[t].vir_lj, vir_lj);
145 void do_force_lowlevel(t_forcerec *fr, t_inputrec *ir,
146 t_idef *idef, t_commrec *cr,
147 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
149 rvec x[], history_t *hist,
152 gmx_enerdata_t *enerd,
173 real dvdl_dum[efptNR], dvdl_nb[efptNR];
176 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
179 set_pbc(&pbc, fr->ePBC, box);
181 /* reset free energy components */
182 for (i = 0; i < efptNR; i++)
189 for (i = 0; (i < DIM); i++)
191 box_size[i] = box[i][i];
196 /* do QMMM first if requested */
199 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
202 /* Call the short range functions all in one go. */
205 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
206 #define TAKETIME FALSE
209 MPI_Barrier(cr->mpi_comm_mygroup);
216 /* foreign lambda component for walls */
217 real dvdl_walls = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
218 enerd->grpp.ener[egLJSR], nrnb);
219 enerd->dvdl_lin[efptVDW] += dvdl_walls;
222 /* If doing GB, reset dvda and calculate the Born radii */
223 if (ir->implicit_solvent)
225 wallcycle_sub_start(wcycle, ewcsNONBONDED);
227 for (i = 0; i < born->nr; i++)
234 calc_gb_rad(cr, fr, ir, top, x, &(fr->gblist), born, md, nrnb);
237 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
241 /* We only do non-bonded calculation with group scheme here, the verlet
242 * calls are done from do_force_cutsVERLET(). */
243 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
246 /* Add short-range interactions */
247 donb_flags |= GMX_NONBONDED_DO_SR;
249 /* Currently all group scheme kernels always calculate (shift-)forces */
250 if (flags & GMX_FORCE_FORCES)
252 donb_flags |= GMX_NONBONDED_DO_FORCE;
254 if (flags & GMX_FORCE_VIRIAL)
256 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
258 if (flags & GMX_FORCE_ENERGY)
260 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
262 if (flags & GMX_FORCE_DO_LR)
264 donb_flags |= GMX_NONBONDED_DO_LR;
267 wallcycle_sub_start(wcycle, ewcsNONBONDED);
268 do_nonbonded(fr, x, f, f_longrange, md, excl,
270 lambda, dvdl_nb, -1, -1, donb_flags);
272 /* If we do foreign lambda and we have soft-core interactions
273 * we have to recalculate the (non-linear) energies contributions.
275 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
277 for (i = 0; i < enerd->n_lambda; i++)
281 for (j = 0; j < efptNR; j++)
283 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
285 reset_foreign_enerdata(enerd);
286 do_nonbonded(fr, x, f, f_longrange, md, excl,
287 &(enerd->foreign_grpp), nrnb,
288 lam_i, dvdl_dum, -1, -1,
289 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
290 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
291 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
294 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
298 /* If we are doing GB, calculate bonded forces and apply corrections
299 * to the solvation forces */
300 /* MRS: Eventually, many need to include free energy contribution here! */
301 if (ir->implicit_solvent)
303 wallcycle_sub_start(wcycle, ewcsLISTED);
304 calc_gb_forces(cr, md, born, top, x, f, fr, idef,
305 ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
306 wallcycle_sub_stop(wcycle, ewcsLISTED);
317 if (fepvals->sc_alpha != 0)
319 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
323 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
326 if (fepvals->sc_alpha != 0)
328 /* even though coulomb part is linear, we already added it, beacuse we
329 need to go through the vdw calculation anyway */
331 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
335 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
343 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
346 /* Shift the coordinates. Must be done before listed forces and PPPM,
347 * but is also necessary for SHAKE and update, therefore it can NOT
348 * go when no listed forces have to be evaluated.
351 /* Here sometimes we would not need to shift with NBFonly,
352 * but we do so anyhow for consistency of the returned coordinates.
356 shift_self(graph, box, x);
359 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
363 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
366 /* Check whether we need to do listed interactions or correct for exclusions */
368 ((flags & GMX_FORCE_LISTED)
369 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype)))
371 /* Since all atoms are in the rectangular or triclinic unit-cell,
372 * only single box vector shifts (2 in x) are required.
374 set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
378 do_force_listed(wcycle, box, ir->fepvals, cr->ms,
379 idef, (const rvec *) x, hist, f, fr,
380 &pbc, graph, enerd, nrnb, lambda, md, fcd,
381 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL,
387 clear_mat(fr->vir_el_recip);
388 clear_mat(fr->vir_lj_recip);
390 /* Do long-range electrostatics and/or LJ-PME, including related short-range
393 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
396 real Vlr_q = 0, Vlr_lj = 0, Vcorr_q = 0, Vcorr_lj = 0;
397 real dvdl_long_range_q = 0, dvdl_long_range_lj = 0;
399 bSB = (ir->nwall == 2);
403 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
404 box_size[ZZ] *= ir->wall_ewald_zfac;
407 if (EEL_PME_EWALD(fr->eeltype) || EVDW_PME(fr->vdwtype))
409 real dvdl_long_range_correction_q = 0;
410 real dvdl_long_range_correction_lj = 0;
411 /* With the Verlet scheme exclusion forces are calculated
412 * in the non-bonded kernel.
414 /* The TPI molecule does not have exclusions with the rest
415 * of the system and no intra-molecular PME grid
416 * contributions will be calculated in
417 * gmx_pme_calc_energy.
419 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
420 ir->ewald_geometry != eewg3D ||
421 ir->epsilon_surface != 0)
425 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
429 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
432 nthreads = gmx_omp_nthreads_get(emntBonded);
433 #pragma omp parallel for num_threads(nthreads) schedule(static)
434 for (t = 0; t < nthreads; t++)
438 tensor *vir_q, *vir_lj;
439 real *Vcorrt_q, *Vcorrt_lj, *dvdlt_q, *dvdlt_lj;
442 fnv = fr->f_novirsum;
443 vir_q = &fr->vir_el_recip;
444 vir_lj = &fr->vir_lj_recip;
446 Vcorrt_lj = &Vcorr_lj;
447 dvdlt_q = &dvdl_long_range_correction_q;
448 dvdlt_lj = &dvdl_long_range_correction_lj;
453 vir_q = &fr->f_t[t].vir_q;
454 vir_lj = &fr->f_t[t].vir_lj;
455 Vcorrt_q = &fr->f_t[t].Vcorr_q;
456 Vcorrt_lj = &fr->f_t[t].Vcorr_lj;
457 dvdlt_q = &fr->f_t[t].dvdl[efptCOUL];
458 dvdlt_lj = &fr->f_t[t].dvdl[efptVDW];
459 for (i = 0; i < fr->natoms_force; i++)
469 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
472 md->nChargePerturbed ? md->chargeB : NULL,
474 md->nTypePerturbed ? md->sqrt_c6B : NULL,
476 md->nTypePerturbed ? md->sigmaB : NULL,
478 md->nTypePerturbed ? md->sigma3B : NULL,
479 ir->cutoff_scheme != ecutsVERLET,
480 excl, x, bSB ? boxs : box, mu_tot,
483 fnv, *vir_q, *vir_lj,
485 lambda[efptCOUL], lambda[efptVDW],
490 reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
491 fr->vir_el_recip, fr->vir_lj_recip,
493 &dvdl_long_range_correction_q,
494 &dvdl_long_range_correction_lj,
497 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
500 if (EEL_PME_EWALD(fr->eeltype) && fr->n_tpi == 0)
502 Vcorr_q += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
503 &dvdl_long_range_correction_q,
507 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction_q;
508 enerd->dvdl_lin[efptVDW] += dvdl_long_range_correction_lj;
510 if ((EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype)) && (cr->duty & DUTY_PME))
512 /* Do reciprocal PME for Coulomb and/or LJ. */
513 assert(fr->n_tpi >= 0);
514 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
516 pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
517 if (EEL_PME(fr->eeltype))
519 pme_flags |= GMX_PME_DO_COULOMB;
521 if (EVDW_PME(fr->vdwtype))
523 pme_flags |= GMX_PME_DO_LJ;
525 if (flags & GMX_FORCE_FORCES)
527 pme_flags |= GMX_PME_CALC_F;
529 if (flags & GMX_FORCE_VIRIAL)
531 pme_flags |= GMX_PME_CALC_ENER_VIR;
535 /* We don't calculate f, but we do want the potential */
536 pme_flags |= GMX_PME_CALC_POT;
538 wallcycle_start(wcycle, ewcPMEMESH);
539 status = gmx_pme_do(fr->pmedata,
540 0, md->homenr - fr->n_tpi,
542 md->chargeA, md->chargeB,
543 md->sqrt_c6A, md->sqrt_c6B,
544 md->sigmaA, md->sigmaB,
545 bSB ? boxs : box, cr,
546 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
547 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
549 fr->vir_el_recip, fr->ewaldcoeff_q,
550 fr->vir_lj_recip, fr->ewaldcoeff_lj,
552 lambda[efptCOUL], lambda[efptVDW],
553 &dvdl_long_range_q, &dvdl_long_range_lj, pme_flags);
554 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
557 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
559 /* We should try to do as little computation after
560 * this as possible, because parallel PME synchronizes
561 * the nodes, so we want all load imbalance of the
562 * rest of the force calculation to be before the PME
563 * call. DD load balancing is done on the whole time
564 * of the force call (without PME).
569 if (EVDW_PME(ir->vdwtype))
572 gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
574 /* Determine the PME grid energy of the test molecule
575 * with the PME grid potential of the other charges.
577 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
578 x + md->homenr - fr->n_tpi,
579 md->chargeA + md->homenr - fr->n_tpi,
585 if (!EEL_PME(fr->eeltype) && EEL_PME_EWALD(fr->eeltype))
587 Vlr_q = do_ewald(ir, x, fr->f_novirsum,
588 md->chargeA, md->chargeB,
589 box_size, cr, md->homenr,
590 fr->vir_el_recip, fr->ewaldcoeff_q,
591 lambda[efptCOUL], &dvdl_long_range_q, fr->ewald_table);
594 /* Note that with separate PME nodes we get the real energies later */
595 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_q;
596 enerd->dvdl_lin[efptVDW] += dvdl_long_range_lj;
597 enerd->term[F_COUL_RECIP] = Vlr_q + Vcorr_q;
598 enerd->term[F_LJ_RECIP] = Vlr_lj + Vcorr_lj;
601 fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
602 Vlr_q, Vcorr_q, enerd->term[F_COUL_RECIP]);
603 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
604 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
605 fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
606 Vlr_lj, Vcorr_lj, enerd->term[F_LJ_RECIP]);
607 pr_rvecs(debug, 0, "vir_lj_recip after corr", fr->vir_lj_recip, DIM);
612 /* Is there a reaction-field exclusion correction needed? */
613 if (EEL_RF(fr->eeltype) && eelRF_NEC != fr->eeltype)
615 /* With the Verlet scheme, exclusion forces are calculated
616 * in the non-bonded kernel.
618 if (ir->cutoff_scheme != ecutsVERLET)
620 real dvdl_rf_excl = 0;
621 enerd->term[F_RF_EXCL] =
622 RF_excl_correction(fr, graph, md, excl, x, f,
623 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
625 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
634 print_nrnb(debug, nrnb);
642 MPI_Barrier(cr->mpi_comm_mygroup);
645 if (fr->timesteps == 11)
648 fprintf(stderr, "* PP load balancing info: rank %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
649 cr->nodeid, gmx_step_str(fr->timesteps, buf),
650 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
651 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
659 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
664 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
668 for (i = 0; i < F_NRE; i++)
671 enerd->foreign_term[i] = 0;
675 for (i = 0; i < efptNR; i++)
677 enerd->dvdl_lin[i] = 0;
678 enerd->dvdl_nonlin[i] = 0;
684 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
686 enerd->grpp.nener = n2;
687 enerd->foreign_grpp.nener = n2;
688 for (i = 0; (i < egNR); i++)
690 snew(enerd->grpp.ener[i], n2);
691 snew(enerd->foreign_grpp.ener[i], n2);
696 enerd->n_lambda = 1 + n_lambda;
697 snew(enerd->enerpart_lambda, enerd->n_lambda);
705 void destroy_enerdata(gmx_enerdata_t *enerd)
709 for (i = 0; (i < egNR); i++)
711 sfree(enerd->grpp.ener[i]);
714 for (i = 0; (i < egNR); i++)
716 sfree(enerd->foreign_grpp.ener[i]);
721 sfree(enerd->enerpart_lambda);
725 static real sum_v(int n, real v[])
731 for (i = 0; (i < n); i++)
739 void sum_epot(gmx_grppairener_t *grpp, real *epot)
743 /* Accumulate energies */
744 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
745 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
746 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
747 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
748 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
749 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
750 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
751 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
753 /* lattice part of LR doesnt belong to any group
754 * and has been added earlier
756 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
757 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
760 for (i = 0; (i < F_EPOT); i++)
762 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
764 epot[F_EPOT] += epot[i];
769 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
774 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
775 enerd->term[F_DVDL] = 0.0;
776 for (i = 0; i < efptNR; i++)
778 if (fepvals->separate_dvdl[i])
780 /* could this be done more readably/compactly? */
793 index = F_DVDL_BONDED;
795 case (efptRESTRAINT):
796 index = F_DVDL_RESTRAINT;
802 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
805 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
806 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
811 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
814 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
815 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
820 /* Notes on the foreign lambda free energy difference evaluation:
821 * Adding the potential and ekin terms that depend linearly on lambda
822 * as delta lam * dvdl to the energy differences is exact.
823 * For the constraints this is not exact, but we have no other option
824 * without literally changing the lengths and reevaluating the energies at each step.
825 * (try to remedy this post 4.6 - MRS)
826 * For the non-bonded LR term we assume that the soft-core (if present)
827 * no longer affects the energy beyond the short-range cut-off,
828 * which is a very good approximation (except for exotic settings).
829 * (investigate how to overcome this post 4.6 - MRS)
831 if (fepvals->separate_dvdl[efptBONDED])
833 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
837 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
839 enerd->term[F_DVDL_CONSTR] = 0;
841 for (i = 0; i < fepvals->n_lambda; i++)
843 /* note we are iterating over fepvals here!
844 For the current lam, dlam = 0 automatically,
845 so we don't need to add anything to the
846 enerd->enerpart_lambda[0] */
848 /* we don't need to worry about dvdl_lin contributions to dE at
849 current lambda, because the contributions to the current
850 lambda are automatically zeroed */
852 for (j = 0; j < efptNR; j++)
854 /* Note that this loop is over all dhdl components, not just the separated ones */
855 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
856 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
859 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
860 fepvals->all_lambda[j][i], efpt_names[j],
861 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
862 dlam, enerd->dvdl_lin[j]);
869 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
873 /* First reset all foreign energy components. Foreign energies always called on
874 neighbor search steps */
875 for (i = 0; (i < egNR); i++)
877 for (j = 0; (j < enerd->grpp.nener); j++)
879 enerd->foreign_grpp.ener[i][j] = 0.0;
883 /* potential energy components */
884 for (i = 0; (i <= F_EPOT); i++)
886 enerd->foreign_term[i] = 0.0;
890 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
891 gmx_enerdata_t *enerd,
897 /* First reset all energy components, except for the long range terms
898 * on the master at non neighbor search steps, since the long range
899 * terms have already been summed at the last neighbor search step.
901 bKeepLR = (fr->bTwinRange && !bNS);
902 for (i = 0; (i < egNR); i++)
904 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
906 for (j = 0; (j < enerd->grpp.nener); j++)
908 enerd->grpp.ener[i][j] = 0.0;
912 for (i = 0; i < efptNR; i++)
914 enerd->dvdl_lin[i] = 0.0;
915 enerd->dvdl_nonlin[i] = 0.0;
918 /* Normal potential energy components */
919 for (i = 0; (i <= F_EPOT); i++)
921 enerd->term[i] = 0.0;
923 /* Initialize the dVdlambda term with the long range contribution */
924 /* Initialize the dvdl term with the long range contribution */
925 enerd->term[F_DVDL] = 0.0;
926 enerd->term[F_DVDL_COUL] = 0.0;
927 enerd->term[F_DVDL_VDW] = 0.0;
928 enerd->term[F_DVDL_BONDED] = 0.0;
929 enerd->term[F_DVDL_RESTRAINT] = 0.0;
930 enerd->term[F_DKDL] = 0.0;
931 if (enerd->n_lambda > 0)
933 for (i = 0; i < enerd->n_lambda; i++)
935 enerd->enerpart_lambda[i] = 0.0;
938 /* reset foreign energy data - separate function since we also call it elsewhere */
939 reset_foreign_enerdata(enerd);