Moved domdec structures out of commrec.h.
[gromacs.git] / src / gromacs / mdlib / force.cpp
blob34d6473503b1f2e58f6b7d5e6f196fbdef7be38e
1 /*
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,2015, 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.
37 #include "gmxpre.h"
39 #include "force.h"
41 #include "config.h"
43 #include <assert.h>
44 #include <math.h>
45 #include <string.h>
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/ewald/ewald.h"
50 #include "gromacs/ewald/long-range-correction.h"
51 #include "gromacs/ewald/pme.h"
52 #include "gromacs/fileio/txtdump.h"
53 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
54 #include "gromacs/gmxlib/network.h"
55 #include "gromacs/gmxlib/nrnb.h"
56 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
57 #include "gromacs/legacyheaders/types/commrec.h"
58 #include "gromacs/listed-forces/listed-forces.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdlib/forcerec-threading.h"
61 #include "gromacs/mdlib/genborn.h"
62 #include "gromacs/mdlib/mdrun.h"
63 #include "gromacs/mdlib/ns.h"
64 #include "gromacs/mdlib/qmmm.h"
65 #include "gromacs/mdtypes/md_enums.h"
66 #include "gromacs/pbcutil/ishift.h"
67 #include "gromacs/pbcutil/mshift.h"
68 #include "gromacs/pbcutil/pbc.h"
69 #include "gromacs/timing/wallcycle.h"
70 #include "gromacs/utility/cstringutil.h"
71 #include "gromacs/utility/exceptions.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/smalloc.h"
75 void ns(FILE *fp,
76 t_forcerec *fr,
77 matrix box,
78 gmx_groups_t *groups,
79 gmx_localtop_t *top,
80 t_mdatoms *md,
81 t_commrec *cr,
82 t_nrnb *nrnb,
83 gmx_bool bFillGrid,
84 gmx_bool bDoLongRangeNS)
86 int nsearch;
89 if (!fr->ns->nblist_initialized)
91 init_neighbor_list(fp, fr, md->homenr);
94 if (fr->bTwinRange)
96 fr->nlr = 0;
99 nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
100 bFillGrid, bDoLongRangeNS);
101 if (debug)
103 fprintf(debug, "nsearch = %d\n", nsearch);
106 /* Check whether we have to do dynamic load balancing */
107 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
108 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
109 &(top->idef),opts->ngener);
111 if (fr->ns->dump_nl > 0)
113 dump_nblist(fp, cr, fr, fr->ns->dump_nl);
117 static void reduce_thread_energies(tensor vir_q, tensor vir_lj,
118 real *Vcorr_q, real *Vcorr_lj,
119 real *dvdl_q, real *dvdl_lj,
120 int nthreads,
121 ewald_corr_thread_t *ewc_t)
123 int t;
125 for (t = 1; t < nthreads; t++)
127 *Vcorr_q += ewc_t[t].Vcorr_q;
128 *Vcorr_lj += ewc_t[t].Vcorr_lj;
129 *dvdl_q += ewc_t[t].dvdl[efptCOUL];
130 *dvdl_lj += ewc_t[t].dvdl[efptVDW];
131 m_add(vir_q, ewc_t[t].vir_q, vir_q);
132 m_add(vir_lj, ewc_t[t].vir_lj, vir_lj);
136 void do_force_lowlevel(t_forcerec *fr, t_inputrec *ir,
137 t_idef *idef, t_commrec *cr,
138 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
139 t_mdatoms *md,
140 rvec x[], history_t *hist,
141 rvec f[],
142 rvec f_longrange[],
143 gmx_enerdata_t *enerd,
144 t_fcdata *fcd,
145 gmx_localtop_t *top,
146 gmx_genborn_t *born,
147 gmx_bool bBornRadii,
148 matrix box,
149 t_lambda *fepvals,
150 real *lambda,
151 t_graph *graph,
152 t_blocka *excl,
153 rvec mu_tot[],
154 int flags,
155 float *cycles_pme)
157 int i, j;
158 int donb_flags;
159 gmx_bool bSB;
160 int pme_flags;
161 matrix boxs;
162 rvec box_size;
163 t_pbc pbc;
164 real dvdl_dum[efptNR], dvdl_nb[efptNR];
166 #ifdef GMX_MPI
167 double t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
168 #endif
170 set_pbc(&pbc, fr->ePBC, box);
172 /* reset free energy components */
173 for (i = 0; i < efptNR; i++)
175 dvdl_nb[i] = 0;
176 dvdl_dum[i] = 0;
179 /* Reset box */
180 for (i = 0; (i < DIM); i++)
182 box_size[i] = box[i][i];
185 /* do QMMM first if requested */
186 if (fr->bQMMM)
188 enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
191 /* Call the short range functions all in one go. */
193 #ifdef GMX_MPI
194 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
195 #define TAKETIME FALSE
196 if (TAKETIME)
198 MPI_Barrier(cr->mpi_comm_mygroup);
199 t0 = MPI_Wtime();
201 #endif
203 if (ir->nwall)
205 /* foreign lambda component for walls */
206 real dvdl_walls = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
207 enerd->grpp.ener[egLJSR], nrnb);
208 enerd->dvdl_lin[efptVDW] += dvdl_walls;
211 /* If doing GB, reset dvda and calculate the Born radii */
212 if (ir->implicit_solvent)
214 wallcycle_sub_start(wcycle, ewcsNONBONDED);
216 for (i = 0; i < born->nr; i++)
218 fr->dvda[i] = 0;
221 if (bBornRadii)
223 calc_gb_rad(cr, fr, ir, top, x, fr->gblist, born, md, nrnb);
226 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
229 where();
230 /* We only do non-bonded calculation with group scheme here, the verlet
231 * calls are done from do_force_cutsVERLET(). */
232 if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
234 donb_flags = 0;
235 /* Add short-range interactions */
236 donb_flags |= GMX_NONBONDED_DO_SR;
238 /* Currently all group scheme kernels always calculate (shift-)forces */
239 if (flags & GMX_FORCE_FORCES)
241 donb_flags |= GMX_NONBONDED_DO_FORCE;
243 if (flags & GMX_FORCE_VIRIAL)
245 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
247 if (flags & GMX_FORCE_ENERGY)
249 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
251 if (flags & GMX_FORCE_DO_LR)
253 donb_flags |= GMX_NONBONDED_DO_LR;
256 wallcycle_sub_start(wcycle, ewcsNONBONDED);
257 do_nonbonded(fr, x, f, f_longrange, md, excl,
258 &enerd->grpp, nrnb,
259 lambda, dvdl_nb, -1, -1, donb_flags);
261 /* If we do foreign lambda and we have soft-core interactions
262 * we have to recalculate the (non-linear) energies contributions.
264 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
266 for (i = 0; i < enerd->n_lambda; i++)
268 real lam_i[efptNR];
270 for (j = 0; j < efptNR; j++)
272 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
274 reset_foreign_enerdata(enerd);
275 do_nonbonded(fr, x, f, f_longrange, md, excl,
276 &(enerd->foreign_grpp), nrnb,
277 lam_i, dvdl_dum, -1, -1,
278 (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
279 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
280 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
283 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
284 where();
287 /* If we are doing GB, calculate bonded forces and apply corrections
288 * to the solvation forces */
289 /* MRS: Eventually, many need to include free energy contribution here! */
290 if (ir->implicit_solvent)
292 wallcycle_sub_start(wcycle, ewcsLISTED);
293 calc_gb_forces(cr, md, born, top, x, f, fr, idef,
294 ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
295 wallcycle_sub_stop(wcycle, ewcsLISTED);
298 #ifdef GMX_MPI
299 if (TAKETIME)
301 t1 = MPI_Wtime();
302 fr->t_fnbf += t1-t0;
304 #endif
306 if (fepvals->sc_alpha != 0)
308 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
310 else
312 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
315 if (fepvals->sc_alpha != 0)
317 /* even though coulomb part is linear, we already added it, beacuse we
318 need to go through the vdw calculation anyway */
320 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
322 else
324 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
327 if (debug)
329 pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
332 /* Shift the coordinates. Must be done before listed forces and PPPM,
333 * but is also necessary for SHAKE and update, therefore it can NOT
334 * go when no listed forces have to be evaluated.
336 * The shifting and PBC code is deliberately not timed, since with
337 * the Verlet scheme it only takes non-zero time with triclinic
338 * boxes, and even then the time is around a factor of 100 less
339 * than the next smallest counter.
343 /* Here sometimes we would not need to shift with NBFonly,
344 * but we do so anyhow for consistency of the returned coordinates.
346 if (graph)
348 shift_self(graph, box, x);
349 if (TRICLINIC(box))
351 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
353 else
355 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
358 /* Check whether we need to do listed interactions or correct for exclusions */
359 if (fr->bMolPBC &&
360 ((flags & GMX_FORCE_LISTED)
361 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype)))
363 /* TODO There are no electrostatics methods that require this
364 transformation, when using the Verlet scheme, so update the
365 above conditional. */
366 /* Since all atoms are in the rectangular or triclinic unit-cell,
367 * only single box vector shifts (2 in x) are required.
369 set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
370 TRUE, box);
373 do_force_listed(wcycle, box, ir->fepvals, cr->ms,
374 idef, (const rvec *) x, hist, f, fr,
375 &pbc, graph, enerd, nrnb, lambda, md, fcd,
376 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL,
377 flags);
379 where();
381 *cycles_pme = 0;
382 clear_mat(fr->vir_el_recip);
383 clear_mat(fr->vir_lj_recip);
385 /* Do long-range electrostatics and/or LJ-PME, including related short-range
386 * corrections.
388 if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
390 int status = 0;
391 real Vlr_q = 0, Vlr_lj = 0, Vcorr_q = 0, Vcorr_lj = 0;
392 real dvdl_long_range_q = 0, dvdl_long_range_lj = 0;
394 bSB = (ir->nwall == 2);
395 if (bSB)
397 copy_mat(box, boxs);
398 svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
399 box_size[ZZ] *= ir->wall_ewald_zfac;
402 if (EEL_PME_EWALD(fr->eeltype) || EVDW_PME(fr->vdwtype))
404 real dvdl_long_range_correction_q = 0;
405 real dvdl_long_range_correction_lj = 0;
406 /* With the Verlet scheme exclusion forces are calculated
407 * in the non-bonded kernel.
409 /* The TPI molecule does not have exclusions with the rest
410 * of the system and no intra-molecular PME grid
411 * contributions will be calculated in
412 * gmx_pme_calc_energy.
414 if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
415 ir->ewald_geometry != eewg3D ||
416 ir->epsilon_surface != 0)
418 int nthreads, t;
420 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
422 if (fr->n_tpi > 0)
424 gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
427 nthreads = fr->nthread_ewc;
428 #pragma omp parallel for num_threads(nthreads) schedule(static)
429 for (t = 0; t < nthreads; t++)
433 tensor *vir_q, *vir_lj;
434 real *Vcorrt_q, *Vcorrt_lj, *dvdlt_q, *dvdlt_lj;
435 if (t == 0)
437 vir_q = &fr->vir_el_recip;
438 vir_lj = &fr->vir_lj_recip;
439 Vcorrt_q = &Vcorr_q;
440 Vcorrt_lj = &Vcorr_lj;
441 dvdlt_q = &dvdl_long_range_correction_q;
442 dvdlt_lj = &dvdl_long_range_correction_lj;
444 else
446 vir_q = &fr->ewc_t[t].vir_q;
447 vir_lj = &fr->ewc_t[t].vir_lj;
448 Vcorrt_q = &fr->ewc_t[t].Vcorr_q;
449 Vcorrt_lj = &fr->ewc_t[t].Vcorr_lj;
450 dvdlt_q = &fr->ewc_t[t].dvdl[efptCOUL];
451 dvdlt_lj = &fr->ewc_t[t].dvdl[efptVDW];
452 clear_mat(*vir_q);
453 clear_mat(*vir_lj);
455 *dvdlt_q = 0;
456 *dvdlt_lj = 0;
458 /* Threading is only supported with the Verlet cut-off
459 * scheme and then only single particle forces (no
460 * exclusion forces) are calculated, so we can store
461 * the forces in the normal, single fr->f_novirsum array.
463 ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
464 cr, t, fr,
465 md->chargeA, md->chargeB,
466 md->sqrt_c6A, md->sqrt_c6B,
467 md->sigmaA, md->sigmaB,
468 md->sigma3A, md->sigma3B,
469 md->nChargePerturbed || md->nTypePerturbed,
470 ir->cutoff_scheme != ecutsVERLET,
471 excl, x, bSB ? boxs : box, mu_tot,
472 ir->ewald_geometry,
473 ir->epsilon_surface,
474 fr->f_novirsum, *vir_q, *vir_lj,
475 Vcorrt_q, Vcorrt_lj,
476 lambda[efptCOUL], lambda[efptVDW],
477 dvdlt_q, dvdlt_lj);
479 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
481 if (nthreads > 1)
483 reduce_thread_energies(fr->vir_el_recip, fr->vir_lj_recip,
484 &Vcorr_q, &Vcorr_lj,
485 &dvdl_long_range_correction_q,
486 &dvdl_long_range_correction_lj,
487 nthreads, fr->ewc_t);
489 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
492 if (EEL_PME_EWALD(fr->eeltype) && fr->n_tpi == 0)
494 /* This is not in a subcounter because it takes a
495 negligible and constant-sized amount of time */
496 Vcorr_q += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
497 &dvdl_long_range_correction_q,
498 fr->vir_el_recip);
501 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction_q;
502 enerd->dvdl_lin[efptVDW] += dvdl_long_range_correction_lj;
504 if ((EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype)) && (cr->duty & DUTY_PME))
506 /* Do reciprocal PME for Coulomb and/or LJ. */
507 assert(fr->n_tpi >= 0);
508 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
510 pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
511 if (EEL_PME(fr->eeltype))
513 pme_flags |= GMX_PME_DO_COULOMB;
515 if (EVDW_PME(fr->vdwtype))
517 pme_flags |= GMX_PME_DO_LJ;
519 if (flags & GMX_FORCE_FORCES)
521 pme_flags |= GMX_PME_CALC_F;
523 if (flags & GMX_FORCE_VIRIAL)
525 pme_flags |= GMX_PME_CALC_ENER_VIR;
527 if (fr->n_tpi > 0)
529 /* We don't calculate f, but we do want the potential */
530 pme_flags |= GMX_PME_CALC_POT;
532 wallcycle_start(wcycle, ewcPMEMESH);
533 status = gmx_pme_do(fr->pmedata,
534 0, md->homenr - fr->n_tpi,
535 x, fr->f_novirsum,
536 md->chargeA, md->chargeB,
537 md->sqrt_c6A, md->sqrt_c6B,
538 md->sigmaA, md->sigmaB,
539 bSB ? boxs : box, cr,
540 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
541 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
542 nrnb, wcycle,
543 fr->vir_el_recip, fr->ewaldcoeff_q,
544 fr->vir_lj_recip, fr->ewaldcoeff_lj,
545 &Vlr_q, &Vlr_lj,
546 lambda[efptCOUL], lambda[efptVDW],
547 &dvdl_long_range_q, &dvdl_long_range_lj, pme_flags);
548 *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
549 if (status != 0)
551 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
553 /* We should try to do as little computation after
554 * this as possible, because parallel PME synchronizes
555 * the nodes, so we want all load imbalance of the
556 * rest of the force calculation to be before the PME
557 * call. DD load balancing is done on the whole time
558 * of the force call (without PME).
561 if (fr->n_tpi > 0)
563 if (EVDW_PME(ir->vdwtype))
566 gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
568 /* Determine the PME grid energy of the test molecule
569 * with the PME grid potential of the other charges.
571 gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
572 x + md->homenr - fr->n_tpi,
573 md->chargeA + md->homenr - fr->n_tpi,
574 &Vlr_q);
579 if (!EEL_PME(fr->eeltype) && EEL_PME_EWALD(fr->eeltype))
581 Vlr_q = do_ewald(ir, x, fr->f_novirsum,
582 md->chargeA, md->chargeB,
583 box_size, cr, md->homenr,
584 fr->vir_el_recip, fr->ewaldcoeff_q,
585 lambda[efptCOUL], &dvdl_long_range_q, fr->ewald_table);
588 /* Note that with separate PME nodes we get the real energies later */
589 enerd->dvdl_lin[efptCOUL] += dvdl_long_range_q;
590 enerd->dvdl_lin[efptVDW] += dvdl_long_range_lj;
591 enerd->term[F_COUL_RECIP] = Vlr_q + Vcorr_q;
592 enerd->term[F_LJ_RECIP] = Vlr_lj + Vcorr_lj;
593 if (debug)
595 fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
596 Vlr_q, Vcorr_q, enerd->term[F_COUL_RECIP]);
597 pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
598 pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
599 fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
600 Vlr_lj, Vcorr_lj, enerd->term[F_LJ_RECIP]);
601 pr_rvecs(debug, 0, "vir_lj_recip after corr", fr->vir_lj_recip, DIM);
604 else
606 /* Is there a reaction-field exclusion correction needed?
607 * With the Verlet scheme, exclusion forces are calculated
608 * in the non-bonded kernel.
610 if (ir->cutoff_scheme != ecutsVERLET && EEL_RF(fr->eeltype))
612 real dvdl_rf_excl = 0;
613 enerd->term[F_RF_EXCL] =
614 RF_excl_correction(fr, graph, md, excl, x, f,
615 fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
617 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
620 where();
622 if (debug)
624 print_nrnb(debug, nrnb);
627 #ifdef GMX_MPI
628 if (TAKETIME)
630 t2 = MPI_Wtime();
631 MPI_Barrier(cr->mpi_comm_mygroup);
632 t3 = MPI_Wtime();
633 fr->t_wait += t3-t2;
634 if (fr->timesteps == 11)
636 char buf[22];
637 fprintf(stderr, "* PP load balancing info: rank %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
638 cr->nodeid, gmx_step_str(fr->timesteps, buf),
639 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
640 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
642 fr->timesteps++;
644 #endif
646 if (debug)
648 pr_rvecs(debug, 0, "fshift after bondeds", fr->fshift, SHIFTS);
653 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
655 int i, n2;
657 for (i = 0; i < F_NRE; i++)
659 enerd->term[i] = 0;
660 enerd->foreign_term[i] = 0;
664 for (i = 0; i < efptNR; i++)
666 enerd->dvdl_lin[i] = 0;
667 enerd->dvdl_nonlin[i] = 0;
670 n2 = ngener*ngener;
671 if (debug)
673 fprintf(debug, "Creating %d sized group matrix for energies\n", n2);
675 enerd->grpp.nener = n2;
676 enerd->foreign_grpp.nener = n2;
677 for (i = 0; (i < egNR); i++)
679 snew(enerd->grpp.ener[i], n2);
680 snew(enerd->foreign_grpp.ener[i], n2);
683 if (n_lambda)
685 enerd->n_lambda = 1 + n_lambda;
686 snew(enerd->enerpart_lambda, enerd->n_lambda);
688 else
690 enerd->n_lambda = 0;
694 void destroy_enerdata(gmx_enerdata_t *enerd)
696 int i;
698 for (i = 0; (i < egNR); i++)
700 sfree(enerd->grpp.ener[i]);
703 for (i = 0; (i < egNR); i++)
705 sfree(enerd->foreign_grpp.ener[i]);
708 if (enerd->n_lambda)
710 sfree(enerd->enerpart_lambda);
714 static real sum_v(int n, real v[])
716 real t;
717 int i;
719 t = 0.0;
720 for (i = 0; (i < n); i++)
722 t = t + v[i];
725 return t;
728 void sum_epot(gmx_grppairener_t *grpp, real *epot)
730 int i;
732 /* Accumulate energies */
733 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
734 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
735 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
736 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
737 epot[F_COUL_LR] = sum_v(grpp->nener, grpp->ener[egCOULLR]);
738 epot[F_LJ_LR] = sum_v(grpp->nener, grpp->ener[egLJLR]);
739 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
740 epot[F_GBPOL] += sum_v(grpp->nener, grpp->ener[egGB]);
742 /* lattice part of LR doesnt belong to any group
743 * and has been added earlier
745 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
746 epot[F_BHAM_LR] = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
748 epot[F_EPOT] = 0;
749 for (i = 0; (i < F_EPOT); i++)
751 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
753 epot[F_EPOT] += epot[i];
758 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
760 int i, j, index;
761 double dlam;
763 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
764 enerd->term[F_DVDL] = 0.0;
765 for (i = 0; i < efptNR; i++)
767 if (fepvals->separate_dvdl[i])
769 /* could this be done more readably/compactly? */
770 switch (i)
772 case (efptMASS):
773 index = F_DKDL;
774 break;
775 case (efptCOUL):
776 index = F_DVDL_COUL;
777 break;
778 case (efptVDW):
779 index = F_DVDL_VDW;
780 break;
781 case (efptBONDED):
782 index = F_DVDL_BONDED;
783 break;
784 case (efptRESTRAINT):
785 index = F_DVDL_RESTRAINT;
786 break;
787 default:
788 index = F_DVDL;
789 break;
791 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
792 if (debug)
794 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
795 efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
798 else
800 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
801 if (debug)
803 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
804 efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
809 /* Notes on the foreign lambda free energy difference evaluation:
810 * Adding the potential and ekin terms that depend linearly on lambda
811 * as delta lam * dvdl to the energy differences is exact.
812 * For the constraints this is not exact, but we have no other option
813 * without literally changing the lengths and reevaluating the energies at each step.
814 * (try to remedy this post 4.6 - MRS)
815 * For the non-bonded LR term we assume that the soft-core (if present)
816 * no longer affects the energy beyond the short-range cut-off,
817 * which is a very good approximation (except for exotic settings).
818 * (investigate how to overcome this post 4.6 - MRS)
820 if (fepvals->separate_dvdl[efptBONDED])
822 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
824 else
826 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
828 enerd->term[F_DVDL_CONSTR] = 0;
830 for (i = 0; i < fepvals->n_lambda; i++)
832 /* note we are iterating over fepvals here!
833 For the current lam, dlam = 0 automatically,
834 so we don't need to add anything to the
835 enerd->enerpart_lambda[0] */
837 /* we don't need to worry about dvdl_lin contributions to dE at
838 current lambda, because the contributions to the current
839 lambda are automatically zeroed */
841 for (j = 0; j < efptNR; j++)
843 /* Note that this loop is over all dhdl components, not just the separated ones */
844 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
845 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
846 if (debug)
848 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
849 fepvals->all_lambda[j][i], efpt_names[j],
850 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
851 dlam, enerd->dvdl_lin[j]);
858 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
860 int i, j;
862 /* First reset all foreign energy components. Foreign energies always called on
863 neighbor search steps */
864 for (i = 0; (i < egNR); i++)
866 for (j = 0; (j < enerd->grpp.nener); j++)
868 enerd->foreign_grpp.ener[i][j] = 0.0;
872 /* potential energy components */
873 for (i = 0; (i <= F_EPOT); i++)
875 enerd->foreign_term[i] = 0.0;
879 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
880 gmx_enerdata_t *enerd,
881 gmx_bool bMaster)
883 gmx_bool bKeepLR;
884 int i, j;
886 /* First reset all energy components, except for the long range terms
887 * on the master at non neighbor search steps, since the long range
888 * terms have already been summed at the last neighbor search step.
890 bKeepLR = (fr->bTwinRange && !bNS);
891 for (i = 0; (i < egNR); i++)
893 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
895 for (j = 0; (j < enerd->grpp.nener); j++)
897 enerd->grpp.ener[i][j] = 0.0;
901 for (i = 0; i < efptNR; i++)
903 enerd->dvdl_lin[i] = 0.0;
904 enerd->dvdl_nonlin[i] = 0.0;
907 /* Normal potential energy components */
908 for (i = 0; (i <= F_EPOT); i++)
910 enerd->term[i] = 0.0;
912 /* Initialize the dVdlambda term with the long range contribution */
913 /* Initialize the dvdl term with the long range contribution */
914 enerd->term[F_DVDL] = 0.0;
915 enerd->term[F_DVDL_COUL] = 0.0;
916 enerd->term[F_DVDL_VDW] = 0.0;
917 enerd->term[F_DVDL_BONDED] = 0.0;
918 enerd->term[F_DVDL_RESTRAINT] = 0.0;
919 enerd->term[F_DKDL] = 0.0;
920 if (enerd->n_lambda > 0)
922 for (i = 0; i < enerd->n_lambda; i++)
924 enerd->enerpart_lambda[i] = 0.0;
927 /* reset foreign energy data - separate function since we also call it elsewhere */
928 reset_foreign_enerdata(enerd);