5bb62cd5b2c378c10c3e407b1d64eae32eba6c8b
[alexxy/gromacs.git] / src / gromacs / mdlib / force.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "force.h"
41
42 #include <cassert>
43 #include <cmath>
44 #include <cstring>
45
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"
72
73 using gmx::ArrayRef;
74 using gmx::RVec;
75
76 static void clearEwaldThreadOutput(ewald_corr_thread_t* ewc_t)
77 {
78     ewc_t->Vcorr_q                                        = 0;
79     ewc_t->Vcorr_lj                                       = 0;
80     ewc_t->dvdl[FreeEnergyPerturbationCouplingType::Coul] = 0;
81     ewc_t->dvdl[FreeEnergyPerturbationCouplingType::Vdw]  = 0;
82     clear_mat(ewc_t->vir_q);
83     clear_mat(ewc_t->vir_lj);
84 }
85
86 static void reduceEwaldThreadOuput(int nthreads, ewald_corr_thread_t* ewc_t)
87 {
88     ewald_corr_thread_t& dest = ewc_t[0];
89
90     for (int t = 1; t < nthreads; t++)
91     {
92         dest.Vcorr_q += ewc_t[t].Vcorr_q;
93         dest.Vcorr_lj += ewc_t[t].Vcorr_lj;
94         dest.dvdl[FreeEnergyPerturbationCouplingType::Coul] +=
95                 ewc_t[t].dvdl[FreeEnergyPerturbationCouplingType::Coul];
96         dest.dvdl[FreeEnergyPerturbationCouplingType::Vdw] +=
97                 ewc_t[t].dvdl[FreeEnergyPerturbationCouplingType::Vdw];
98         m_add(dest.vir_q, ewc_t[t].vir_q, dest.vir_q);
99         m_add(dest.vir_lj, ewc_t[t].vir_lj, dest.vir_lj);
100     }
101 }
102
103 void calculateLongRangeNonbondeds(t_forcerec*                   fr,
104                                   const t_inputrec&             ir,
105                                   const t_commrec*              cr,
106                                   t_nrnb*                       nrnb,
107                                   gmx_wallcycle_t               wcycle,
108                                   const t_mdatoms*              md,
109                                   gmx::ArrayRef<const RVec>     coordinates,
110                                   gmx::ForceWithVirial*         forceWithVirial,
111                                   gmx_enerdata_t*               enerd,
112                                   const matrix                  box,
113                                   const real*                   lambda,
114                                   const rvec*                   mu_tot,
115                                   const gmx::StepWorkload&      stepWork,
116                                   const DDBalanceRegionHandler& ddBalanceRegionHandler)
117 {
118     const bool computePmeOnCpu = (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
119                                  && thisRankHasDuty(cr, DUTY_PME)
120                                  && (pme_run_mode(fr->pmedata) == PmeRunMode::CPU);
121
122     const bool haveEwaldSurfaceTerm = haveEwaldSurfaceContribution(ir);
123
124     /* Do long-range electrostatics and/or LJ-PME
125      * and compute PME surface terms when necessary.
126      */
127     if ((computePmeOnCpu || fr->ic->eeltype == CoulombInteractionType::Ewald || haveEwaldSurfaceTerm)
128         && stepWork.computeNonbondedForces)
129     {
130         int  status = 0;
131         real Vlr_q = 0, Vlr_lj = 0;
132
133         /* We reduce all virial, dV/dlambda and energy contributions, except
134          * for the reciprocal energies (Vlr_q, Vlr_lj) into the same struct.
135          */
136         ewald_corr_thread_t& ewaldOutput = fr->ewc_t[0];
137         clearEwaldThreadOutput(&ewaldOutput);
138
139         if (EEL_PME_EWALD(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
140         {
141             /* Calculate the Ewald surface force and energy contributions, when necessary */
142             if (haveEwaldSurfaceTerm)
143             {
144                 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
145
146                 int nthreads = fr->nthread_ewc;
147 #pragma omp parallel for num_threads(nthreads) schedule(static)
148                 for (int t = 0; t < nthreads; t++)
149                 {
150                     try
151                     {
152                         ewald_corr_thread_t& ewc_t = fr->ewc_t[t];
153                         if (t > 0)
154                         {
155                             clearEwaldThreadOutput(&ewc_t);
156                         }
157
158                         /* Threading is only supported with the Verlet cut-off
159                          * scheme and then only single particle forces (no
160                          * exclusion forces) are calculated, so we can store
161                          * the forces in the normal, single forceWithVirial->force_ array.
162                          */
163                         const rvec* x = as_rvec_array(coordinates.data());
164                         ewald_LRcorrection(
165                                 md->homenr,
166                                 cr,
167                                 nthreads,
168                                 t,
169                                 *fr,
170                                 ir,
171                                 md->chargeA,
172                                 md->chargeB,
173                                 (md->nChargePerturbed != 0),
174                                 x,
175                                 box,
176                                 mu_tot,
177                                 as_rvec_array(forceWithVirial->force_.data()),
178                                 &ewc_t.Vcorr_q,
179                                 lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
180                                 &ewc_t.dvdl[FreeEnergyPerturbationCouplingType::Coul]);
181                     }
182                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
183                 }
184                 if (nthreads > 1)
185                 {
186                     reduceEwaldThreadOuput(nthreads, fr->ewc_t);
187                 }
188                 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
189             }
190
191             if (EEL_PME_EWALD(fr->ic->eeltype) && fr->n_tpi == 0)
192             {
193                 /* This is not in a subcounter because it takes a
194                    negligible and constant-sized amount of time */
195                 ewaldOutput.Vcorr_q += ewald_charge_correction(
196                         cr,
197                         fr,
198                         lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
199                         box,
200                         &ewaldOutput.dvdl[FreeEnergyPerturbationCouplingType::Coul],
201                         ewaldOutput.vir_q);
202             }
203
204             if (computePmeOnCpu)
205             {
206                 /* Do reciprocal PME for Coulomb and/or LJ. */
207                 assert(fr->n_tpi >= 0);
208                 if (fr->n_tpi == 0 || stepWork.stateChanged)
209                 {
210                     /* With domain decomposition we close the CPU side load
211                      * balancing region here, because PME does global
212                      * communication that acts as a global barrier.
213                      */
214                     ddBalanceRegionHandler.closeAfterForceComputationCpu();
215
216                     wallcycle_start(wcycle, ewcPMEMESH);
217                     status = gmx_pme_do(
218                             fr->pmedata,
219                             gmx::constArrayRefFromArray(coordinates.data(), md->homenr - fr->n_tpi),
220                             forceWithVirial->force_,
221                             md->chargeA,
222                             md->chargeB,
223                             md->sqrt_c6A,
224                             md->sqrt_c6B,
225                             md->sigmaA,
226                             md->sigmaB,
227                             box,
228                             cr,
229                             DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
230                             DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
231                             nrnb,
232                             wcycle,
233                             ewaldOutput.vir_q,
234                             ewaldOutput.vir_lj,
235                             &Vlr_q,
236                             &Vlr_lj,
237                             lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
238                             lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Vdw)],
239                             &ewaldOutput.dvdl[FreeEnergyPerturbationCouplingType::Coul],
240                             &ewaldOutput.dvdl[FreeEnergyPerturbationCouplingType::Vdw],
241                             stepWork);
242                     wallcycle_stop(wcycle, ewcPMEMESH);
243                     if (status != 0)
244                     {
245                         gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
246                     }
247
248                     /* We should try to do as little computation after
249                      * this as possible, because parallel PME synchronizes
250                      * the nodes, so we want all load imbalance of the
251                      * rest of the force calculation to be before the PME
252                      * call.  DD load balancing is done on the whole time
253                      * of the force call (without PME).
254                      */
255                 }
256                 if (fr->n_tpi > 0)
257                 {
258                     /* Determine the PME grid energy of the test molecule
259                      * with the PME grid potential of the other charges.
260                      */
261                     gmx_pme_calc_energy(
262                             fr->pmedata,
263                             coordinates.subArray(md->homenr - fr->n_tpi, fr->n_tpi),
264                             gmx::arrayRefFromArray(md->chargeA + md->homenr - fr->n_tpi, fr->n_tpi),
265                             &Vlr_q);
266                 }
267             }
268         }
269
270         if (fr->ic->eeltype == CoulombInteractionType::Ewald)
271         {
272             const rvec* x = as_rvec_array(coordinates.data());
273             Vlr_q         = do_ewald(ir,
274                              x,
275                              as_rvec_array(forceWithVirial->force_.data()),
276                              md->chargeA,
277                              md->chargeB,
278                              box,
279                              cr,
280                              md->homenr,
281                              ewaldOutput.vir_q,
282                              fr->ic->ewaldcoeff_q,
283                              lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Coul)],
284                              &ewaldOutput.dvdl[FreeEnergyPerturbationCouplingType::Coul],
285                              fr->ewald_table.get());
286         }
287
288         /* Note that with separate PME nodes we get the real energies later */
289         // TODO it would be simpler if we just accumulated a single
290         // long-range virial contribution.
291         forceWithVirial->addVirialContribution(ewaldOutput.vir_q);
292         forceWithVirial->addVirialContribution(ewaldOutput.vir_lj);
293         enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Coul] +=
294                 ewaldOutput.dvdl[FreeEnergyPerturbationCouplingType::Coul];
295         enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Vdw] +=
296                 ewaldOutput.dvdl[FreeEnergyPerturbationCouplingType::Vdw];
297         enerd->term[F_COUL_RECIP] = Vlr_q + ewaldOutput.Vcorr_q;
298         enerd->term[F_LJ_RECIP]   = Vlr_lj + ewaldOutput.Vcorr_lj;
299
300         if (debug)
301         {
302             fprintf(debug,
303                     "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
304                     Vlr_q,
305                     ewaldOutput.Vcorr_q,
306                     enerd->term[F_COUL_RECIP]);
307             pr_rvecs(debug, 0, "vir_el_recip after corr", ewaldOutput.vir_q, DIM);
308             fprintf(debug,
309                     "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
310                     Vlr_lj,
311                     ewaldOutput.Vcorr_lj,
312                     enerd->term[F_LJ_RECIP]);
313             pr_rvecs(debug, 0, "vir_lj_recip after corr", ewaldOutput.vir_lj, DIM);
314         }
315     }
316
317     if (debug)
318     {
319         print_nrnb(debug, nrnb);
320     }
321 }