Simplify vsite PBC handling
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nonbonded.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,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 "nonbonded.h"
40
41 #include <cassert>
42 #include <cstdio>
43 #include <cstdlib>
44
45 #include "thread_mpi/threads.h"
46
47 #include "gromacs/gmxlib/nrnb.h"
48 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
49 #include "gromacs/gmxlib/nonbonded/nb_generic.h"
50 #include "gromacs/gmxlib/nonbonded/nb_generic_cg.h"
51 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
52 #include "gromacs/listed_forces/bonded.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdtypes/enerdata.h"
56 #include "gromacs/mdtypes/forcerec.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/mdtypes/mdatom.h"
59 #include "gromacs/mdtypes/nblist.h"
60 #include "gromacs/pbcutil/ishift.h"
61 #include "gromacs/pbcutil/mshift.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/tables/forcetable.h"
64 #include "gromacs/utility/arraysize.h"
65 #include "gromacs/utility/basedefinitions.h"
66 #include "gromacs/utility/cstringutil.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/smalloc.h"
69
70 void
71 gmx_nonbonded_set_kernel_pointers(t_nblist *nl)
72 {
73     const char *     elec;
74     const char *     elec_mod;
75     const char *     vdw;
76     const char *     vdw_mod;
77     const char *     geom;
78
79     nl->kernelptr_vf = nullptr;
80     nl->kernelptr_v  = nullptr;
81     nl->kernelptr_f  = nullptr;
82
83     elec     = gmx_nbkernel_elec_names[nl->ielec];
84     elec_mod = eintmod_names[nl->ielecmod];
85     vdw      = gmx_nbkernel_vdw_names[nl->ivdw];
86     vdw_mod  = eintmod_names[nl->ivdwmod];
87     geom     = gmx_nblist_geometry_names[nl->igeometry];
88
89     if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
90     {
91         nl->kernelptr_vf       = reinterpret_cast<void *>(gmx_nb_free_energy_kernel);
92         nl->kernelptr_f        = reinterpret_cast<void *>(gmx_nb_free_energy_kernel);
93         nl->simd_padding_width = 1;
94     }
95     else if (!gmx_strcasecmp_min(geom, "CG-CG"))
96     {
97         nl->kernelptr_vf       = reinterpret_cast<void *>(gmx_nb_generic_cg_kernel);
98         nl->kernelptr_f        = reinterpret_cast<void *>(gmx_nb_generic_cg_kernel);
99         nl->simd_padding_width = 1;
100     }
101     else
102     {
103         /* "Pick" the only remaining kernel, the generic one.
104          * We only do this for particle-particle kernels; by leaving the water-optimized kernel
105          * pointers to NULL, the water optimization will automatically be disabled for this interaction.
106          */
107         if (!gmx_strcasecmp_min(geom, "Particle-Particle"))
108         {
109             nl->kernelptr_vf       = reinterpret_cast<void *>(gmx_nb_generic_kernel);
110             nl->kernelptr_f        = reinterpret_cast<void *>(gmx_nb_generic_kernel);
111             nl->simd_padding_width = 1;
112             if (debug)
113             {
114                 fprintf(debug,
115                         "WARNING - Slow generic NB kernel used for neighborlist with\n"
116                         "    Elec: '%s', Modifier: '%s'\n"
117                         "    Vdw:  '%s', Modifier: '%s'\n",
118                         elec, elec_mod, vdw, vdw_mod);
119             }
120         }
121     }
122 }
123
124 void do_nonbonded(const t_forcerec  *fr,
125                   rvec               x[],
126                   rvec               f_shortrange[],
127                   const t_mdatoms   *mdatoms,
128                   const t_blocka    *excl,
129                   gmx_grppairener_t *grppener,
130                   t_nrnb            *nrnb,
131                   real              *lambda,
132                   real              *dvdl,
133                   int                nls,
134                   int                eNL,
135                   int                flags)
136 {
137     t_nblist *        nlist;
138     int               n, n0, n1, i, i0, i1;
139     t_nblists *       nblists;
140     nb_kernel_data_t  kernel_data;
141     nb_kernel_t *     kernelptr = nullptr;
142     rvec *            f;
143
144     kernel_data.flags                   = flags;
145     kernel_data.exclusions              = excl;
146     kernel_data.lambda                  = lambda;
147     kernel_data.dvdl                    = dvdl;
148
149     if (eNL >= 0)
150     {
151         i0 = eNL;
152         i1 = i0+1;
153     }
154     else
155     {
156         i0 = 0;
157         i1 = eNL_NR;
158     }
159
160     if (nls >= 0)
161     {
162         n0 = nls;
163         n1 = nls+1;
164     }
165     else
166     {
167         n0 = 0;
168         n1 = fr->nnblists;
169     }
170
171     for (n = n0; (n < n1); n++)
172     {
173         nblists = &fr->nblists[n];
174
175         /* Tabulated kernels hard-code a lot of assumptions about the
176          * structure of these tables, but that's not worth fixing with
177          * the group scheme due for removal soon. As a token
178          * improvement, this assertion will stop code segfaulting if
179          * someone assumes that extending the group-scheme table-type
180          * enumeration is something that GROMACS supports. */
181         static_assert(etiNR == 3, "");
182
183         kernel_data.table_elec              = nblists->table_elec;
184         kernel_data.table_vdw               = nblists->table_vdw;
185         kernel_data.table_elec_vdw          = nblists->table_elec_vdw;
186
187         {
188             {
189                 /* Short-range */
190                 if (!(flags & GMX_NONBONDED_DO_SR))
191                 {
192                     continue;
193                 }
194                 kernel_data.energygrp_elec          = grppener->ener[egCOULSR];
195                 kernel_data.energygrp_vdw           = grppener->ener[fr->bBHAM ? egBHAMSR : egLJSR];
196                 nlist = nblists->nlist_sr;
197                 f                                   = f_shortrange;
198             }
199
200             for (i = i0; (i < i1); i++)
201             {
202                 if (nlist[i].nri > 0)
203                 {
204                     if (flags & GMX_NONBONDED_DO_POTENTIAL)
205                     {
206                         /* Potential and force */
207                         kernelptr = reinterpret_cast<nb_kernel_t *>(nlist[i].kernelptr_vf);
208                     }
209                     else
210                     {
211                         /* Force only, no potential */
212                         kernelptr = reinterpret_cast<nb_kernel_t *>(nlist[i].kernelptr_f);
213                     }
214
215                     if (nlist[i].type != GMX_NBLIST_INTERACTION_FREE_ENERGY && (flags & GMX_NONBONDED_DO_FOREIGNLAMBDA))
216                     {
217                         /* We don't need the non-perturbed interactions */
218                         continue;
219                     }
220                     /* Neighborlists whose kernelptr==NULL will always be empty */
221                     if (kernelptr != nullptr)
222                     {
223                         (*kernelptr)(&(nlist[i]), x, f, const_cast<t_forcerec*>(fr),
224                                      const_cast<t_mdatoms*>(mdatoms), &kernel_data, nrnb);
225                     }
226                     else
227                     {
228                         gmx_fatal(FARGS, "Non-empty neighborlist does not have any kernel pointer assigned.");
229                     }
230                 }
231             }
232         }
233     }
234 }