2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * \brief This file defines high-level functions for mdrun to compute
38 * energies and forces for listed interactions.
40 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_listed-forces
46 #include "listed-forces.h"
54 #include "gromacs/legacyheaders/disre.h"
55 #include "gromacs/legacyheaders/force.h"
56 #include "gromacs/legacyheaders/network.h"
57 #include "gromacs/legacyheaders/nrnb.h"
58 #include "gromacs/legacyheaders/orires.h"
59 #include "gromacs/legacyheaders/types/force_flags.h"
60 #include "gromacs/listed-forces/bonded.h"
61 #include "gromacs/listed-forces/position-restraints.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/forcerec-threading.h"
64 #include "gromacs/pbcutil/ishift.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/timing/wallcycle.h"
67 #include "gromacs/utility/smalloc.h"
74 /*! \brief Return true if ftype is an explicit pair-listed LJ or
75 * COULOMB interaction type: bonded LJ (usually 1-4), or special
76 * listed non-bonded for FEP. */
78 isPairInteraction(int ftype)
80 return ((ftype) >= F_LJ14 && (ftype) <= F_LJC_PAIRS_NB);
83 /*! \brief Zero thread-local force-output buffers */
85 zero_thread_forces(f_thread_t *f_t, int n,
86 int nblock, int blocksize)
88 int b, a0, a1, a, i, j;
90 if (n > f_t->f_nalloc)
92 f_t->f_nalloc = over_alloc_large(n);
93 srenew(f_t->f, f_t->f_nalloc);
96 if (!bitmask_is_zero(f_t->red_mask))
98 for (b = 0; b < nblock; b++)
100 if (bitmask_is_set(f_t->red_mask, b))
103 a1 = std::min((b+1)*blocksize, n);
104 for (a = a0; a < a1; a++)
106 clear_rvec(f_t->f[a]);
111 for (i = 0; i < SHIFTS; i++)
113 clear_rvec(f_t->fshift[i]);
115 for (i = 0; i < F_NRE; i++)
119 for (i = 0; i < egNR; i++)
121 for (j = 0; j < f_t->grpp.nener; j++)
123 f_t->grpp.ener[i][j] = 0;
126 for (i = 0; i < efptNR; i++)
132 /*! \brief The max thread number is arbitrary, we used a fixed number
133 * to avoid memory management. Using more than 16 threads is probably
134 * never useful performance wise. */
135 #define MAX_BONDED_THREADS 256
137 /*! \brief Reduce thread-local force buffers */
139 reduce_thread_force_buffer(int n, rvec *f,
140 int nthreads, f_thread_t *f_t,
141 int nblock, int block_size)
145 if (nthreads > MAX_BONDED_THREADS)
147 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
151 /* This reduction can run on any number of threads,
152 * independently of nthreads.
154 #pragma omp parallel for num_threads(nthreads) schedule(static)
155 for (b = 0; b < nblock; b++)
157 rvec *fp[MAX_BONDED_THREADS];
161 /* Determine which threads contribute to this block */
163 for (ft = 1; ft < nthreads; ft++)
165 if (bitmask_is_set(f_t[ft].red_mask, b))
167 fp[nfb++] = f_t[ft].f;
172 /* Reduce force buffers for threads that contribute */
174 a1 = (b+1)*block_size;
175 a1 = std::min(a1, n);
176 for (a = a0; a < a1; a++)
178 for (fb = 0; fb < nfb; fb++)
180 rvec_inc(f[a], fp[fb][a]);
187 /*! \brief Reduce thread-local forces */
189 reduce_thread_forces(int n, rvec *f, rvec *fshift,
190 real *ener, gmx_grppairener_t *grpp, real *dvdl,
191 int nthreads, f_thread_t *f_t,
192 int nblock, int block_size,
193 gmx_bool bCalcEnerVir,
198 /* Reduce the bonded force buffer */
199 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
202 /* When necessary, reduce energy and virial using one thread only */
207 for (i = 0; i < SHIFTS; i++)
209 for (t = 1; t < nthreads; t++)
211 rvec_inc(fshift[i], f_t[t].fshift[i]);
214 for (i = 0; i < F_NRE; i++)
216 for (t = 1; t < nthreads; t++)
218 ener[i] += f_t[t].ener[i];
221 for (i = 0; i < egNR; i++)
223 for (j = 0; j < f_t[1].grpp.nener; j++)
225 for (t = 1; t < nthreads; t++)
228 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
234 for (i = 0; i < efptNR; i++)
237 for (t = 1; t < nthreads; t++)
239 dvdl[i] += f_t[t].dvdl[i];
246 /*! \brief Calculate one element of the list of bonded interactions
249 calc_one_bond(int thread,
250 int ftype, const t_idef *idef,
251 const rvec x[], rvec f[], rvec fshift[],
253 const t_pbc *pbc, const t_graph *g,
254 gmx_grppairener_t *grpp,
256 real *lambda, real *dvdl,
257 const t_mdatoms *md, t_fcdata *fcd,
258 gmx_bool bCalcEnerVir,
259 int *global_atom_index)
261 int nat1, nbonds, efptFTYPE;
266 if (IS_RESTRAINT_TYPE(ftype))
268 efptFTYPE = efptRESTRAINT;
272 efptFTYPE = efptBONDED;
275 nat1 = interaction_function[ftype].nratoms + 1;
276 nbonds = idef->il[ftype].nr/nat1;
277 iatoms = idef->il[ftype].iatoms;
279 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
280 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
282 if (!isPairInteraction(ftype))
286 v = cmap_dihs(nbn, iatoms+nb0,
287 idef->iparams, &idef->cmap_grid,
289 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
290 md, fcd, global_atom_index);
292 #ifdef GMX_SIMD_HAVE_REAL
293 else if (ftype == F_ANGLES &&
294 !bCalcEnerVir && fr->efep == efepNO)
296 /* No energies, shift forces, dvdl */
297 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
300 pbc, g, lambda[efptFTYPE], md, fcd,
305 else if (ftype == F_PDIHS &&
306 !bCalcEnerVir && fr->efep == efepNO)
308 /* No energies, shift forces, dvdl */
309 #ifdef GMX_SIMD_HAVE_REAL
314 (nbn, idef->il[ftype].iatoms+nb0,
317 pbc, g, lambda[efptFTYPE], md, fcd,
323 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
326 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
327 md, fcd, global_atom_index);
332 v = do_pairs(ftype, nbn, iatoms+nb0, idef->iparams, x, f, fshift,
333 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
338 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
346 void calc_listed(const gmx_multisim_t *ms,
348 const rvec x[], history_t *hist,
349 rvec f[], t_forcerec *fr,
350 const struct t_pbc *pbc,
351 const struct t_pbc *pbc_full,
352 const struct t_graph *g,
353 gmx_enerdata_t *enerd, t_nrnb *nrnb,
356 t_fcdata *fcd, int *global_atom_index,
359 gmx_bool bCalcEnerVir;
361 real dvdl[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
362 of lambda, which will be thrown away in the end*/
363 const t_pbc *pbc_null;
366 assert(fr->nthreads == idef->nthreads);
368 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
370 for (i = 0; i < efptNR; i++)
386 p_graph(debug, "Bondage is fun", g);
390 if (idef->il[F_POSRES].nr > 0)
392 posres_wrapper(nrnb, idef, pbc_full, x, enerd, lambda, fr);
395 if (idef->il[F_FBPOSRES].nr > 0)
397 fbposres_wrapper(nrnb, idef, pbc_full, x, enerd, fr);
400 /* Do pre force calculation stuff which might require communication */
401 if (idef->il[F_ORIRES].nr)
403 enerd->term[F_ORIRESDEV] =
404 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
405 idef->il[F_ORIRES].iatoms,
406 idef->iparams, md, x,
407 pbc_null, fcd, hist);
409 if (idef->il[F_DISRES].nr)
411 calc_disres_R_6(idef->il[F_DISRES].nr,
412 idef->il[F_DISRES].iatoms,
413 idef->iparams, x, pbc_null,
416 if (fcd->disres.nsystems > 1)
418 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
423 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
424 for (thread = 0; thread < fr->nthreads; thread++)
431 gmx_grppairener_t *grpp;
443 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
444 fr->red_nblock, 1<<fr->red_ashift);
446 ft = fr->f_t[thread].f;
447 fshift = fr->f_t[thread].fshift;
448 epot = fr->f_t[thread].ener;
449 grpp = &fr->f_t[thread].grpp;
450 dvdlt = fr->f_t[thread].dvdl;
452 /* Loop over all bonded force types to calculate the bonded forces */
453 for (ftype = 0; (ftype < F_NRE); ftype++)
455 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
457 v = calc_one_bond(thread, ftype, idef, x,
458 ft, fshift, fr, pbc_null, g, grpp,
460 md, fcd, bCalcEnerVir,
466 if (fr->nthreads > 1)
468 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
469 enerd->term, &enerd->grpp, dvdl,
470 fr->nthreads, fr->f_t,
471 fr->red_nblock, 1<<fr->red_ashift,
473 force_flags & GMX_FORCE_DHDL);
475 if (force_flags & GMX_FORCE_DHDL)
477 for (i = 0; i < efptNR; i++)
479 enerd->dvdl_nonlin[i] += dvdl[i];
483 /* Copy the sum of violations for the distance restraints from fcd */
486 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
491 void calc_listed_lambda(const t_idef *idef,
494 const struct t_pbc *pbc, const struct t_graph *g,
495 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
499 int *global_atom_index)
501 int ftype, nr_nonperturbed, nr;
503 real dvdl_dum[efptNR];
505 const t_pbc *pbc_null;
517 /* Copy the whole idef, so we can modify the contents locally */
519 idef_fe.nthreads = 1;
520 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
522 /* We already have the forces, so we use temp buffers here */
523 snew(f, fr->natoms_force);
524 snew(fshift, SHIFTS);
526 /* Loop over all bonded force types to calculate the bonded energies */
527 for (ftype = 0; (ftype < F_NRE); ftype++)
529 if (ftype_is_bonded_potential(ftype))
531 /* Set the work range of thread 0 to the perturbed bondeds only */
532 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
533 nr = idef->il[ftype].nr;
534 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
535 idef_fe.il_thread_division[ftype*2+1] = nr;
537 /* This is only to get the flop count correct */
538 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
540 if (nr - nr_nonperturbed > 0)
542 v = calc_one_bond(0, ftype, &idef_fe,
543 x, f, fshift, fr, pbc_null, g,
544 grpp, nrnb, lambda, dvdl_dum,
555 sfree(idef_fe.il_thread_division);
559 do_force_listed(gmx_wallcycle *wcycle,
561 const t_lambda *fepvals,
562 const gmx_multisim_t *ms,
568 const struct t_pbc *pbc,
569 const struct t_graph *graph,
570 gmx_enerdata_t *enerd,
575 int *global_atom_index,
578 t_pbc pbc_full; /* Full PBC is needed for position restraints */
580 if (!(flags & GMX_FORCE_LISTED))
585 wallcycle_sub_start(wcycle, ewcsLISTED);
587 if ((idef->il[F_POSRES].nr > 0) ||
588 (idef->il[F_FBPOSRES].nr > 0))
590 set_pbc(&pbc_full, fr->ePBC, box);
592 calc_listed(ms, idef, x, hist, f, fr, pbc, &pbc_full,
593 graph, enerd, nrnb, lambda, md, fcd,
594 global_atom_index, flags);
596 /* Check if we have to determine energy differences
597 * at foreign lambda's.
599 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL))
601 posres_wrapper_lambda(fepvals, idef, &pbc_full, x, enerd, lambda, fr);
603 if (idef->ilsort != ilsortNO_FE)
605 if (idef->ilsort != ilsortFE_SORTED)
607 gmx_incons("The bonded interactions are not sorted for free energy");
609 for (int i = 0; i < enerd->n_lambda; i++)
613 reset_foreign_enerdata(enerd);
614 for (int j = 0; j < efptNR; j++)
616 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
618 calc_listed_lambda(idef, x, fr, pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
619 fcd, global_atom_index);
620 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
621 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
627 wallcycle_sub_stop(wcycle, ewcsLISTED);