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,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "gromacs/domdec/dlbtiming.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/gmxlib/network.h"
53 #include "gromacs/gmxlib/nrnb.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/math/vecdump.h"
56 #include "gromacs/mdlib/forcerec_threading.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/mdtypes/enerdata.h"
59 #include "gromacs/mdtypes/forceoutput.h"
60 #include "gromacs/mdtypes/forcerec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/interaction_const.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/mdtypes/mdatom.h"
65 #include "gromacs/mdtypes/simulation_workload.h"
66 #include "gromacs/pbcutil/ishift.h"
67 #include "gromacs/pbcutil/pbc.h"
68 #include "gromacs/timing/wallcycle.h"
69 #include "gromacs/utility/exceptions.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/smalloc.h"
76 static void clearEwaldThreadOutput(ewald_corr_thread_t* ewc_t)
80 ewc_t->dvdl[efptCOUL] = 0;
81 ewc_t->dvdl[efptVDW] = 0;
82 clear_mat(ewc_t->vir_q);
83 clear_mat(ewc_t->vir_lj);
86 static void reduceEwaldThreadOuput(int nthreads, ewald_corr_thread_t* ewc_t)
88 ewald_corr_thread_t& dest = ewc_t[0];
90 for (int t = 1; t < nthreads; t++)
92 dest.Vcorr_q += ewc_t[t].Vcorr_q;
93 dest.Vcorr_lj += ewc_t[t].Vcorr_lj;
94 dest.dvdl[efptCOUL] += ewc_t[t].dvdl[efptCOUL];
95 dest.dvdl[efptVDW] += ewc_t[t].dvdl[efptVDW];
96 m_add(dest.vir_q, ewc_t[t].vir_q, dest.vir_q);
97 m_add(dest.vir_lj, ewc_t[t].vir_lj, dest.vir_lj);
101 void calculateLongRangeNonbondeds(t_forcerec* fr,
102 const t_inputrec& ir,
105 gmx_wallcycle_t wcycle,
107 gmx::ArrayRef<const RVec> coordinates,
108 gmx::ForceWithVirial* forceWithVirial,
109 gmx_enerdata_t* enerd,
113 const gmx::StepWorkload& stepWork,
114 const DDBalanceRegionHandler& ddBalanceRegionHandler)
116 const bool computePmeOnCpu = (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
117 && thisRankHasDuty(cr, DUTY_PME)
118 && (pme_run_mode(fr->pmedata) == PmeRunMode::CPU);
120 const bool haveEwaldSurfaceTerm = haveEwaldSurfaceContribution(ir);
122 /* Do long-range electrostatics and/or LJ-PME
123 * and compute PME surface terms when necessary.
125 if ((computePmeOnCpu || fr->ic->eeltype == eelEWALD || haveEwaldSurfaceTerm)
126 && stepWork.computeNonbondedForces)
129 real Vlr_q = 0, Vlr_lj = 0;
131 /* We reduce all virial, dV/dlambda and energy contributions, except
132 * for the reciprocal energies (Vlr_q, Vlr_lj) into the same struct.
134 ewald_corr_thread_t& ewaldOutput = fr->ewc_t[0];
135 clearEwaldThreadOutput(&ewaldOutput);
137 if (EEL_PME_EWALD(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
139 /* Calculate the Ewald surface force and energy contributions, when necessary */
140 if (haveEwaldSurfaceTerm)
142 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
144 int nthreads = fr->nthread_ewc;
145 #pragma omp parallel for num_threads(nthreads) schedule(static)
146 for (int t = 0; t < nthreads; t++)
150 ewald_corr_thread_t& ewc_t = fr->ewc_t[t];
153 clearEwaldThreadOutput(&ewc_t);
156 /* Threading is only supported with the Verlet cut-off
157 * scheme and then only single particle forces (no
158 * exclusion forces) are calculated, so we can store
159 * the forces in the normal, single forceWithVirial->force_ array.
161 const rvec* x = as_rvec_array(coordinates.data());
162 ewald_LRcorrection(md->homenr,
170 (md->nChargePerturbed != 0),
174 as_rvec_array(forceWithVirial->force_.data()),
177 &ewc_t.dvdl[efptCOUL]);
179 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
183 reduceEwaldThreadOuput(nthreads, fr->ewc_t);
185 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
188 if (EEL_PME_EWALD(fr->ic->eeltype) && fr->n_tpi == 0)
190 /* This is not in a subcounter because it takes a
191 negligible and constant-sized amount of time */
192 ewaldOutput.Vcorr_q += ewald_charge_correction(
193 cr, fr, lambda[efptCOUL], box, &ewaldOutput.dvdl[efptCOUL], ewaldOutput.vir_q);
198 /* Do reciprocal PME for Coulomb and/or LJ. */
199 assert(fr->n_tpi >= 0);
200 if (fr->n_tpi == 0 || stepWork.stateChanged)
202 /* With domain decomposition we close the CPU side load
203 * balancing region here, because PME does global
204 * communication that acts as a global barrier.
206 ddBalanceRegionHandler.closeAfterForceComputationCpu();
208 wallcycle_start(wcycle, ewcPMEMESH);
211 gmx::constArrayRefFromArray(coordinates.data(), md->homenr - fr->n_tpi),
212 forceWithVirial->force_,
221 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
222 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
231 &ewaldOutput.dvdl[efptCOUL],
232 &ewaldOutput.dvdl[efptVDW],
234 wallcycle_stop(wcycle, ewcPMEMESH);
237 gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
240 /* We should try to do as little computation after
241 * this as possible, because parallel PME synchronizes
242 * the nodes, so we want all load imbalance of the
243 * rest of the force calculation to be before the PME
244 * call. DD load balancing is done on the whole time
245 * of the force call (without PME).
250 /* Determine the PME grid energy of the test molecule
251 * with the PME grid potential of the other charges.
255 coordinates.subArray(md->homenr - fr->n_tpi, fr->n_tpi),
256 gmx::arrayRefFromArray(md->chargeA + md->homenr - fr->n_tpi, fr->n_tpi),
262 if (fr->ic->eeltype == eelEWALD)
264 const rvec* x = as_rvec_array(coordinates.data());
267 as_rvec_array(forceWithVirial->force_.data()),
274 fr->ic->ewaldcoeff_q,
276 &ewaldOutput.dvdl[efptCOUL],
280 /* Note that with separate PME nodes we get the real energies later */
281 // TODO it would be simpler if we just accumulated a single
282 // long-range virial contribution.
283 forceWithVirial->addVirialContribution(ewaldOutput.vir_q);
284 forceWithVirial->addVirialContribution(ewaldOutput.vir_lj);
285 enerd->dvdl_lin[efptCOUL] += ewaldOutput.dvdl[efptCOUL];
286 enerd->dvdl_lin[efptVDW] += ewaldOutput.dvdl[efptVDW];
287 enerd->term[F_COUL_RECIP] = Vlr_q + ewaldOutput.Vcorr_q;
288 enerd->term[F_LJ_RECIP] = Vlr_lj + ewaldOutput.Vcorr_lj;
293 "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
296 enerd->term[F_COUL_RECIP]);
297 pr_rvecs(debug, 0, "vir_el_recip after corr", ewaldOutput.vir_q, DIM);
299 "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
301 ewaldOutput.Vcorr_lj,
302 enerd->term[F_LJ_RECIP]);
303 pr_rvecs(debug, 0, "vir_lj_recip after corr", ewaldOutput.vir_lj, DIM);
309 print_nrnb(debug, nrnb);