3e7c359e0326f478548ef83dc85e7523e3af3729
[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,2018,2019, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "force.h"
40
41 #include <cassert>
42 #include <cmath>
43 #include <cstring>
44
45 #include "gromacs/domdec/dlbtiming.h"
46 #include "gromacs/domdec/domdec.h"
47 #include "gromacs/domdec/domdec_struct.h"
48 #include "gromacs/ewald/ewald.h"
49 #include "gromacs/ewald/long_range_correction.h"
50 #include "gromacs/ewald/pme.h"
51 #include "gromacs/gmxlib/network.h"
52 #include "gromacs/gmxlib/nrnb.h"
53 #include "gromacs/listed_forces/listed_forces.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/math/vecdump.h"
56 #include "gromacs/mdlib/forcerec_threading.h"
57 #include "gromacs/mdlib/qmmm.h"
58 #include "gromacs/mdlib/rf_util.h"
59 #include "gromacs/mdlib/wall.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/enerdata.h"
62 #include "gromacs/mdtypes/forceoutput.h"
63 #include "gromacs/mdtypes/forcerec.h"
64 #include "gromacs/mdtypes/inputrec.h"
65 #include "gromacs/mdtypes/md_enums.h"
66 #include "gromacs/mdtypes/mdatom.h"
67 #include "gromacs/mdtypes/simulation_workload.h"
68 #include "gromacs/pbcutil/ishift.h"
69 #include "gromacs/pbcutil/mshift.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/timing/wallcycle.h"
72 #include "gromacs/utility/exceptions.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/smalloc.h"
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[efptCOUL] = 0;
81     ewc_t->dvdl[efptVDW]  = 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[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);
98     }
99 }
100
101 void do_force_lowlevel(t_forcerec*                         fr,
102                        const t_inputrec*                   ir,
103                        const t_idef*                       idef,
104                        const t_commrec*                    cr,
105                        const gmx_multisim_t*               ms,
106                        t_nrnb*                             nrnb,
107                        gmx_wallcycle_t                     wcycle,
108                        const t_mdatoms*                    md,
109                        gmx::ArrayRefWithPadding<gmx::RVec> coordinates,
110                        history_t*                          hist,
111                        gmx::ForceOutputs*                  forceOutputs,
112                        gmx_enerdata_t*                     enerd,
113                        t_fcdata*                           fcd,
114                        const matrix                        box,
115                        const real*                         lambda,
116                        const t_graph*                      graph,
117                        const rvec*                         mu_tot,
118                        const gmx::StepWorkload&            stepWork,
119                        const DDBalanceRegionHandler&       ddBalanceRegionHandler)
120 {
121     // TODO: Replace all uses of x by const coordinates
122     rvec* x = as_rvec_array(coordinates.paddedArrayRef().data());
123
124     auto& forceWithVirial = forceOutputs->forceWithVirial();
125
126     /* do QMMM first if requested */
127     if (fr->bQMMM)
128     {
129         enerd->term[F_EQM] = calculate_QMMM(cr, &forceOutputs->forceWithShiftForces(), fr->qr);
130     }
131
132     /* Call the short range functions all in one go. */
133
134     if (ir->nwall)
135     {
136         /* foreign lambda component for walls */
137         real dvdl_walls = do_walls(*ir, *fr, box, *md, x, &forceWithVirial, lambda[efptVDW],
138                                    enerd->grpp.ener[egLJSR].data(), nrnb);
139         enerd->dvdl_lin[efptVDW] += dvdl_walls;
140     }
141
142     /* Shift the coordinates. Must be done before listed forces and PPPM,
143      * but is also necessary for SHAKE and update, therefore it can NOT
144      * go when no listed forces have to be evaluated.
145      *
146      * The shifting and PBC code is deliberately not timed, since with
147      * the Verlet scheme it only takes non-zero time with triclinic
148      * boxes, and even then the time is around a factor of 100 less
149      * than the next smallest counter.
150      */
151
152
153     /* Here sometimes we would not need to shift with NBFonly,
154      * but we do so anyhow for consistency of the returned coordinates.
155      */
156     if (graph)
157     {
158         shift_self(graph, box, x);
159         if (TRICLINIC(box))
160         {
161             inc_nrnb(nrnb, eNR_SHIFTX, 2 * graph->nnodes);
162         }
163         else
164         {
165             inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
166         }
167     }
168
169     {
170         t_pbc pbc;
171
172         /* Check whether we need to take into account PBC in listed interactions. */
173         const auto needPbcForListedForces =
174                 fr->bMolPBC && stepWork.computeListedForces && haveCpuListedForces(*fr, *idef, *fcd);
175         if (needPbcForListedForces)
176         {
177             /* Since all atoms are in the rectangular or triclinic unit-cell,
178              * only single box vector shifts (2 in x) are required.
179              */
180             set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd->numCells : nullptr, TRUE, box);
181         }
182
183         do_force_listed(wcycle, box, ir->fepvals, cr, ms, idef, x, hist, forceOutputs, fr, &pbc,
184                         graph, enerd, nrnb, lambda, md, fcd,
185                         DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr, stepWork);
186     }
187
188     const bool computePmeOnCpu = (EEL_PME(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
189                                  && thisRankHasDuty(cr, DUTY_PME)
190                                  && (pme_run_mode(fr->pmedata) == PmeRunMode::CPU);
191
192     const bool haveEwaldSurfaceTerm = haveEwaldSurfaceContribution(*ir);
193
194     /* Do long-range electrostatics and/or LJ-PME
195      * and compute PME surface terms when necessary.
196      */
197     if (computePmeOnCpu || fr->ic->eeltype == eelEWALD || haveEwaldSurfaceTerm)
198     {
199         int  status = 0;
200         real Vlr_q = 0, Vlr_lj = 0;
201
202         /* We reduce all virial, dV/dlambda and energy contributions, except
203          * for the reciprocal energies (Vlr_q, Vlr_lj) into the same struct.
204          */
205         ewald_corr_thread_t& ewaldOutput = fr->ewc_t[0];
206         clearEwaldThreadOutput(&ewaldOutput);
207
208         if (EEL_PME_EWALD(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
209         {
210             /* Calculate the Ewald surface force and energy contributions, when necessary */
211             if (haveEwaldSurfaceTerm)
212             {
213                 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
214
215                 if (fr->n_tpi > 0)
216                 {
217                     gmx_fatal(FARGS,
218                               "TPI with PME currently only works in a 3D geometry with tin-foil "
219                               "boundary conditions");
220                 }
221
222                 int nthreads = fr->nthread_ewc;
223 #pragma omp parallel for num_threads(nthreads) schedule(static)
224                 for (int t = 0; t < nthreads; t++)
225                 {
226                     try
227                     {
228                         ewald_corr_thread_t& ewc_t = fr->ewc_t[t];
229                         if (t > 0)
230                         {
231                             clearEwaldThreadOutput(&ewc_t);
232                         }
233
234                         /* Threading is only supported with the Verlet cut-off
235                          * scheme and then only single particle forces (no
236                          * exclusion forces) are calculated, so we can store
237                          * the forces in the normal, single forceWithVirial->force_ array.
238                          */
239                         ewald_LRcorrection(md->homenr, cr, nthreads, t, *fr, *ir, md->chargeA,
240                                            md->chargeB, (md->nChargePerturbed != 0), x, box, mu_tot,
241                                            as_rvec_array(forceWithVirial.force_.data()),
242                                            &ewc_t.Vcorr_q, lambda[efptCOUL], &ewc_t.dvdl[efptCOUL]);
243                     }
244                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
245                 }
246                 if (nthreads > 1)
247                 {
248                     reduceEwaldThreadOuput(nthreads, fr->ewc_t);
249                 }
250                 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
251             }
252
253             if (EEL_PME_EWALD(fr->ic->eeltype) && fr->n_tpi == 0)
254             {
255                 /* This is not in a subcounter because it takes a
256                    negligible and constant-sized amount of time */
257                 ewaldOutput.Vcorr_q += ewald_charge_correction(
258                         cr, fr, lambda[efptCOUL], box, &ewaldOutput.dvdl[efptCOUL], ewaldOutput.vir_q);
259             }
260
261             if (computePmeOnCpu)
262             {
263                 /* Do reciprocal PME for Coulomb and/or LJ. */
264                 assert(fr->n_tpi >= 0);
265                 if (fr->n_tpi == 0 || stepWork.stateChanged)
266                 {
267                     int pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
268
269                     if (stepWork.computeForces)
270                     {
271                         pme_flags |= GMX_PME_CALC_F;
272                     }
273                     if (stepWork.computeVirial)
274                     {
275                         pme_flags |= GMX_PME_CALC_ENER_VIR;
276                     }
277                     if (fr->n_tpi > 0)
278                     {
279                         /* We don't calculate f, but we do want the potential */
280                         pme_flags |= GMX_PME_CALC_POT;
281                     }
282
283                     /* With domain decomposition we close the CPU side load
284                      * balancing region here, because PME does global
285                      * communication that acts as a global barrier.
286                      */
287                     ddBalanceRegionHandler.closeAfterForceComputationCpu();
288
289                     wallcycle_start(wcycle, ewcPMEMESH);
290                     status = gmx_pme_do(
291                             fr->pmedata,
292                             gmx::constArrayRefFromArray(coordinates.unpaddedConstArrayRef().data(),
293                                                         md->homenr - fr->n_tpi),
294                             forceWithVirial.force_, md->chargeA, md->chargeB, md->sqrt_c6A,
295                             md->sqrt_c6B, md->sigmaA, md->sigmaB, box, cr,
296                             DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
297                             DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0, nrnb, wcycle,
298                             ewaldOutput.vir_q, ewaldOutput.vir_lj, &Vlr_q, &Vlr_lj,
299                             lambda[efptCOUL], lambda[efptVDW], &ewaldOutput.dvdl[efptCOUL],
300                             &ewaldOutput.dvdl[efptVDW], pme_flags);
301                     wallcycle_stop(wcycle, ewcPMEMESH);
302                     if (status != 0)
303                     {
304                         gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
305                     }
306
307                     /* We should try to do as little computation after
308                      * this as possible, because parallel PME synchronizes
309                      * the nodes, so we want all load imbalance of the
310                      * rest of the force calculation to be before the PME
311                      * call.  DD load balancing is done on the whole time
312                      * of the force call (without PME).
313                      */
314                 }
315                 if (fr->n_tpi > 0)
316                 {
317                     if (EVDW_PME(ir->vdwtype))
318                     {
319
320                         gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
321                     }
322                     /* Determine the PME grid energy of the test molecule
323                      * with the PME grid potential of the other charges.
324                      */
325                     gmx_pme_calc_energy(
326                             fr->pmedata,
327                             coordinates.unpaddedConstArrayRef().subArray(md->homenr - fr->n_tpi, fr->n_tpi),
328                             gmx::arrayRefFromArray(md->chargeA + md->homenr - fr->n_tpi, fr->n_tpi),
329                             &Vlr_q);
330                 }
331             }
332         }
333
334         if (fr->ic->eeltype == eelEWALD)
335         {
336             Vlr_q = do_ewald(ir, x, as_rvec_array(forceWithVirial.force_.data()), md->chargeA,
337                              md->chargeB, box, cr, md->homenr, ewaldOutput.vir_q, fr->ic->ewaldcoeff_q,
338                              lambda[efptCOUL], &ewaldOutput.dvdl[efptCOUL], fr->ewald_table);
339         }
340
341         /* Note that with separate PME nodes we get the real energies later */
342         // TODO it would be simpler if we just accumulated a single
343         // long-range virial contribution.
344         forceWithVirial.addVirialContribution(ewaldOutput.vir_q);
345         forceWithVirial.addVirialContribution(ewaldOutput.vir_lj);
346         enerd->dvdl_lin[efptCOUL] += ewaldOutput.dvdl[efptCOUL];
347         enerd->dvdl_lin[efptVDW] += ewaldOutput.dvdl[efptVDW];
348         enerd->term[F_COUL_RECIP] = Vlr_q + ewaldOutput.Vcorr_q;
349         enerd->term[F_LJ_RECIP]   = Vlr_lj + ewaldOutput.Vcorr_lj;
350
351         if (debug)
352         {
353             fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n", Vlr_q,
354                     ewaldOutput.Vcorr_q, enerd->term[F_COUL_RECIP]);
355             pr_rvecs(debug, 0, "vir_el_recip after corr", ewaldOutput.vir_q, DIM);
356             rvec* fshift = as_rvec_array(forceOutputs->forceWithShiftForces().shiftForces().data());
357             pr_rvecs(debug, 0, "fshift after LR Corrections", fshift, SHIFTS);
358             fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n", Vlr_lj,
359                     ewaldOutput.Vcorr_lj, enerd->term[F_LJ_RECIP]);
360             pr_rvecs(debug, 0, "vir_lj_recip after corr", ewaldOutput.vir_lj, DIM);
361         }
362     }
363
364     if (debug)
365     {
366         print_nrnb(debug, nrnb);
367     }
368
369     if (debug)
370     {
371         rvec* fshift = as_rvec_array(forceOutputs->forceWithShiftForces().shiftForces().data());
372         pr_rvecs(debug, 0, "fshift after bondeds", fshift, SHIFTS);
373     }
374 }