Tidy: modernize-use-nullptr
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nonbonded.cpp
blob2d683c34cb9a636229660b3c049383540bc50344
1 /*
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,2017, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "nonbonded.h"
41 #include "config.h"
43 #include <cassert>
44 #include <cstdio>
45 #include <cstdlib>
47 #include "thread_mpi/threads.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
51 #include "gromacs/gmxlib/nonbonded/nb_generic.h"
52 #include "gromacs/gmxlib/nonbonded/nb_generic_cg.h"
53 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
54 #include "gromacs/listed-forces/bonded.h"
55 #include "gromacs/math/utilities.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdtypes/forcerec.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/mdtypes/mdatom.h"
60 #include "gromacs/mdtypes/nblist.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/pbcutil/mshift.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/tables/forcetable.h"
65 #include "gromacs/utility/arraysize.h"
66 #include "gromacs/utility/basedefinitions.h"
67 #include "gromacs/utility/cstringutil.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/smalloc.h"
71 /* Different default (c) and SIMD instructions interaction-specific kernels */
72 #include "gromacs/gmxlib/nonbonded/nb_kernel_c/nb_kernel_c.h"
74 #if GMX_SIMD_X86_SSE2 && !GMX_DOUBLE
75 # include "gromacs/gmxlib/nonbonded/nb_kernel_sse2_single/nb_kernel_sse2_single.h"
76 #endif
77 #if GMX_SIMD_X86_SSE4_1 && !GMX_DOUBLE
78 # include "gromacs/gmxlib/nonbonded/nb_kernel_sse4_1_single/nb_kernel_sse4_1_single.h"
79 #endif
80 #if GMX_SIMD_X86_AVX_128_FMA && !GMX_DOUBLE
81 # include "gromacs/gmxlib/nonbonded/nb_kernel_avx_128_fma_single/nb_kernel_avx_128_fma_single.h"
82 #endif
83 #if (GMX_SIMD_X86_AVX_256 || GMX_SIMD_X86_AVX2_256) && !GMX_DOUBLE
84 # include "gromacs/gmxlib/nonbonded/nb_kernel_avx_256_single/nb_kernel_avx_256_single.h"
85 #endif
86 #if GMX_SIMD_X86_SSE2 && GMX_DOUBLE
87 # include "gromacs/gmxlib/nonbonded/nb_kernel_sse2_double/nb_kernel_sse2_double.h"
88 #endif
89 #if GMX_SIMD_X86_SSE4_1 && GMX_DOUBLE
90 # include "gromacs/gmxlib/nonbonded/nb_kernel_sse4_1_double/nb_kernel_sse4_1_double.h"
91 #endif
92 #if GMX_SIMD_X86_AVX_128_FMA && GMX_DOUBLE
93 # include "gromacs/gmxlib/nonbonded/nb_kernel_avx_128_fma_double/nb_kernel_avx_128_fma_double.h"
94 #endif
95 #if (GMX_SIMD_X86_AVX_256 || GMX_SIMD_X86_AVX2_256) && GMX_DOUBLE
96 # include "gromacs/gmxlib/nonbonded/nb_kernel_avx_256_double/nb_kernel_avx_256_double.h"
97 #endif
98 #if GMX_SIMD_SPARC64_HPC_ACE && GMX_DOUBLE
99 # include "gromacs/gmxlib/nonbonded/nb_kernel_sparc64_hpc_ace_double/nb_kernel_sparc64_hpc_ace_double.h"
100 #endif
103 static tMPI_Thread_mutex_t nonbonded_setup_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
104 static gmx_bool nonbonded_setup_done = FALSE;
107 void
108 gmx_nonbonded_setup(t_forcerec * fr,
109 gmx_bool bGenericKernelOnly)
111 tMPI_Thread_mutex_lock(&nonbonded_setup_mutex);
112 /* Here we are guaranteed only one thread made it. */
113 if (nonbonded_setup_done == FALSE)
115 if (bGenericKernelOnly == FALSE)
117 /* Add the generic kernels to the structure stored statically in nb_kernel.c */
118 nb_kernel_list_add_kernels(kernellist_c, kernellist_c_size);
120 if (!(fr != nullptr && fr->use_simd_kernels == FALSE))
122 /* Add interaction-specific kernels for different architectures */
123 /* Single precision */
124 #if GMX_SIMD_X86_SSE2 && !GMX_DOUBLE
125 nb_kernel_list_add_kernels(kernellist_sse2_single, kernellist_sse2_single_size);
126 #endif
127 #if GMX_SIMD_X86_SSE4_1 && !GMX_DOUBLE
128 nb_kernel_list_add_kernels(kernellist_sse4_1_single, kernellist_sse4_1_single_size);
129 #endif
130 #if GMX_SIMD_X86_AVX_128_FMA && !GMX_DOUBLE
131 nb_kernel_list_add_kernels(kernellist_avx_128_fma_single, kernellist_avx_128_fma_single_size);
132 #endif
133 #if (GMX_SIMD_X86_AVX_256 || GMX_SIMD_X86_AVX2_256) && !GMX_DOUBLE
134 nb_kernel_list_add_kernels(kernellist_avx_256_single, kernellist_avx_256_single_size);
135 #endif
136 /* Double precision */
137 #if GMX_SIMD_X86_SSE2 && GMX_DOUBLE
138 nb_kernel_list_add_kernels(kernellist_sse2_double, kernellist_sse2_double_size);
139 #endif
140 #if GMX_SIMD_X86_SSE4_1 && GMX_DOUBLE
141 nb_kernel_list_add_kernels(kernellist_sse4_1_double, kernellist_sse4_1_double_size);
142 #endif
143 #if GMX_SIMD_X86_AVX_128_FMA && GMX_DOUBLE
144 nb_kernel_list_add_kernels(kernellist_avx_128_fma_double, kernellist_avx_128_fma_double_size);
145 #endif
146 #if (GMX_SIMD_X86_AVX_256 || GMX_SIMD_X86_AVX2_256) && GMX_DOUBLE
147 nb_kernel_list_add_kernels(kernellist_avx_256_double, kernellist_avx_256_double_size);
148 #endif
149 #if GMX_SIMD_SPARC64_HPC_ACE && GMX_DOUBLE
150 nb_kernel_list_add_kernels(kernellist_sparc64_hpc_ace_double, kernellist_sparc64_hpc_ace_double_size);
151 #endif
152 ; /* empty statement to avoid a completely empty block */
155 /* Create a hash for faster lookups */
156 nb_kernel_list_hash_init();
158 nonbonded_setup_done = TRUE;
160 tMPI_Thread_mutex_unlock(&nonbonded_setup_mutex);
165 void
166 gmx_nonbonded_set_kernel_pointers(FILE *log, t_nblist *nl, gmx_bool bElecAndVdwSwitchDiffers)
168 const char * elec;
169 const char * elec_mod;
170 const char * vdw;
171 const char * vdw_mod;
172 const char * geom;
173 const char * other;
175 struct
177 const char * arch;
178 int simd_padding_width;
180 arch_and_padding[] =
182 /* Single precision */
183 #if (GMX_SIMD_X86_AVX_256 || GMX_SIMD_X86_AVX2_256) && !GMX_DOUBLE
184 { "avx_256_single", 8 },
185 #endif
186 #if GMX_SIMD_X86_AVX_128_FMA && !GMX_DOUBLE
187 { "avx_128_fma_single", 4 },
188 #endif
189 #if GMX_SIMD_X86_SSE4_1 && !GMX_DOUBLE
190 { "sse4_1_single", 4 },
191 #endif
192 #if GMX_SIMD_X86_SSE2 && !GMX_DOUBLE
193 { "sse2_single", 4 },
194 #endif
195 /* Double precision */
196 #if (GMX_SIMD_X86_AVX_256 || GMX_SIMD_X86_AVX2_256) && GMX_DOUBLE
197 { "avx_256_double", 4 },
198 #endif
199 #if GMX_SIMD_X86_AVX_128_FMA && GMX_DOUBLE
200 /* Sic. Double precision 2-way SIMD does not require neighbor list padding,
201 * since the kernels execute a loop unrolled a factor 2, followed by
202 * a possible single odd-element epilogue.
204 { "avx_128_fma_double", 1 },
205 #endif
206 #if GMX_SIMD_X86_SSE2 && GMX_DOUBLE
207 /* No padding - see comment above */
208 { "sse2_double", 1 },
209 #endif
210 #if GMX_SIMD_X86_SSE4_1 && GMX_DOUBLE
211 /* No padding - see comment above */
212 { "sse4_1_double", 1 },
213 #endif
214 #if GMX_SIMD_SPARC64_HPC_ACE && GMX_DOUBLE
215 /* No padding - see comment above */
216 { "sparc64_hpc_ace_double", 1 },
217 #endif
218 { "c", 1 },
220 int narch = asize(arch_and_padding);
221 int i;
223 if (nonbonded_setup_done == FALSE)
225 /* We typically call this setup routine before starting timers,
226 * but if that has not been done for whatever reason we do it now.
228 gmx_nonbonded_setup(nullptr, FALSE);
231 /* Not used yet */
232 other = "";
234 nl->kernelptr_vf = nullptr;
235 nl->kernelptr_v = nullptr;
236 nl->kernelptr_f = nullptr;
238 elec = gmx_nbkernel_elec_names[nl->ielec];
239 elec_mod = eintmod_names[nl->ielecmod];
240 vdw = gmx_nbkernel_vdw_names[nl->ivdw];
241 vdw_mod = eintmod_names[nl->ivdwmod];
242 geom = gmx_nblist_geometry_names[nl->igeometry];
244 if (nl->type == GMX_NBLIST_INTERACTION_FREE_ENERGY)
246 nl->kernelptr_vf = (void *) gmx_nb_free_energy_kernel;
247 nl->kernelptr_f = (void *) gmx_nb_free_energy_kernel;
248 nl->simd_padding_width = 1;
250 else if (!gmx_strcasecmp_min(geom, "CG-CG"))
252 nl->kernelptr_vf = (void *) gmx_nb_generic_cg_kernel;
253 nl->kernelptr_f = (void *) gmx_nb_generic_cg_kernel;
254 nl->simd_padding_width = 1;
256 else
258 /* Try to find a specific kernel first */
260 for (i = 0; i < narch && nl->kernelptr_vf == nullptr; i++)
262 nl->kernelptr_vf = (void *) nb_kernel_list_findkernel(log, arch_and_padding[i].arch, elec, elec_mod, vdw, vdw_mod, geom, other, "PotentialAndForce");
263 nl->simd_padding_width = arch_and_padding[i].simd_padding_width;
265 for (i = 0; i < narch && nl->kernelptr_f == nullptr; i++)
267 nl->kernelptr_f = (void *) nb_kernel_list_findkernel(log, arch_and_padding[i].arch, elec, elec_mod, vdw, vdw_mod, geom, other, "Force");
268 nl->simd_padding_width = arch_and_padding[i].simd_padding_width;
270 /* If there is not force-only optimized kernel, is there a potential & force one? */
271 if (nl->kernelptr_f == nullptr)
273 nl->kernelptr_f = (void *) nb_kernel_list_findkernel(nullptr, arch_and_padding[i].arch, elec, elec_mod, vdw, vdw_mod, geom, other, "PotentialAndForce");
274 nl->simd_padding_width = arch_and_padding[i].simd_padding_width;
278 /* For now, the accelerated kernels cannot handle the combination of switch functions for both
279 * electrostatics and VdW that use different switch radius or switch cutoff distances
280 * (both of them enter in the switch function calculation). This would require
281 * us to evaluate two completely separate switch functions for every interaction.
282 * Instead, we disable such kernels by setting the pointer to NULL.
283 * This will cause the generic kernel (which can handle it) to be called instead.
285 * Note that we typically already enable tabulated coulomb interactions for this case,
286 * so this is mostly a safe-guard to make sure we call the generic kernel if the
287 * tables are disabled.
289 if ((nl->ielec != GMX_NBKERNEL_ELEC_NONE) && (nl->ielecmod == eintmodPOTSWITCH) &&
290 (nl->ivdw != GMX_NBKERNEL_VDW_NONE) && (nl->ivdwmod == eintmodPOTSWITCH) &&
291 bElecAndVdwSwitchDiffers)
293 nl->kernelptr_vf = nullptr;
294 nl->kernelptr_f = nullptr;
297 /* Give up, pick a generic one instead.
298 * We only do this for particle-particle kernels; by leaving the water-optimized kernel
299 * pointers to NULL, the water optimization will automatically be disabled for this interaction.
301 if (nl->kernelptr_vf == nullptr && !gmx_strcasecmp_min(geom, "Particle-Particle"))
303 nl->kernelptr_vf = (void *) gmx_nb_generic_kernel;
304 nl->kernelptr_f = (void *) gmx_nb_generic_kernel;
305 nl->simd_padding_width = 1;
306 if (debug)
308 fprintf(debug,
309 "WARNING - Slow generic NB kernel used for neighborlist with\n"
310 " Elec: '%s', Modifier: '%s'\n"
311 " Vdw: '%s', Modifier: '%s'\n",
312 elec, elec_mod, vdw, vdw_mod);
316 return;
319 void do_nonbonded(t_forcerec *fr,
320 rvec x[], rvec f_shortrange[], t_mdatoms *mdatoms, t_blocka *excl,
321 gmx_grppairener_t *grppener,
322 t_nrnb *nrnb, real *lambda, real *dvdl,
323 int nls, int eNL, int flags)
325 t_nblist * nlist;
326 int n, n0, n1, i, i0, i1;
327 t_nblists * nblists;
328 nb_kernel_data_t kernel_data;
329 nb_kernel_t * kernelptr = nullptr;
330 rvec * f;
332 kernel_data.flags = flags;
333 kernel_data.exclusions = excl;
334 kernel_data.lambda = lambda;
335 kernel_data.dvdl = dvdl;
337 if (fr->bAllvsAll)
339 gmx_incons("All-vs-all kernels have not been implemented in version 4.6");
340 return;
343 if (eNL >= 0)
345 i0 = eNL;
346 i1 = i0+1;
348 else
350 i0 = 0;
351 i1 = eNL_NR;
354 if (nls >= 0)
356 n0 = nls;
357 n1 = nls+1;
359 else
361 n0 = 0;
362 n1 = fr->nnblists;
365 for (n = n0; (n < n1); n++)
367 nblists = &fr->nblists[n];
369 /* Tabulated kernels hard-code a lot of assumptions about the
370 * structure of these tables, but that's not worth fixing with
371 * the group scheme due for removal soon. As a token
372 * improvement, this assertion will stop code segfaulting if
373 * someone assumes that extending the group-scheme table-type
374 * enumeration is something that GROMACS supports. */
375 /* cppcheck-suppress duplicateExpression */
376 assert(etiNR == 3);
378 kernel_data.table_elec = nblists->table_elec;
379 kernel_data.table_vdw = nblists->table_vdw;
380 kernel_data.table_elec_vdw = nblists->table_elec_vdw;
384 /* Short-range */
385 if (!(flags & GMX_NONBONDED_DO_SR))
387 continue;
389 kernel_data.energygrp_elec = grppener->ener[egCOULSR];
390 kernel_data.energygrp_vdw = grppener->ener[fr->bBHAM ? egBHAMSR : egLJSR];
391 kernel_data.energygrp_polarization = grppener->ener[egGB];
392 nlist = nblists->nlist_sr;
393 f = f_shortrange;
396 for (i = i0; (i < i1); i++)
398 if (nlist[i].nri > 0)
400 if (flags & GMX_NONBONDED_DO_POTENTIAL)
402 /* Potential and force */
403 kernelptr = (nb_kernel_t *)nlist[i].kernelptr_vf;
405 else
407 /* Force only, no potential */
408 kernelptr = (nb_kernel_t *)nlist[i].kernelptr_f;
411 if (nlist[i].type != GMX_NBLIST_INTERACTION_FREE_ENERGY && (flags & GMX_NONBONDED_DO_FOREIGNLAMBDA))
413 /* We don't need the non-perturbed interactions */
414 continue;
416 /* Neighborlists whose kernelptr==NULL will always be empty */
417 if (kernelptr != nullptr)
419 (*kernelptr)(&(nlist[i]), x, f, fr, mdatoms, &kernel_data, nrnb);
421 else
423 gmx_fatal(FARGS, "Non-empty neighborlist does not have any kernel pointer assigned.");