Fix OpenCL compiles on Mac OS
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_ocl / nbnxn_ocl_jit_support.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2017,2018, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  *  \brief Defines functions that support JIT compilation (e.g. for OpenCL)
37  *
38  *  \author Dimitrios Karkoulis <dimitris.karkoulis@gmail.com>
39  *  \author Mark Abraham <mark.j.abraham@gmail.com>
40  *  \ingroup module_mdlib
41  */
42 #include "gmxpre.h"
43
44 #include <stdlib.h>
45
46 #include <cassert>
47
48 #include <string>
49
50 #include "gromacs/gpu_utils/gpu_utils.h"
51 #include "gromacs/gpu_utils/ocl_compiler.h"
52 #include "gromacs/mdlib/nbnxn_consts.h"
53 #include "gromacs/mdlib/nbnxn_gpu.h"
54 #include "gromacs/mdlib/nbnxn_gpu_jit_support.h"
55 #include "gromacs/mdtypes/interaction_const.h"
56 #include "gromacs/mdtypes/md_enums.h"
57 #include "gromacs/pbcutil/ishift.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/exceptions.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/stringutil.h"
62
63 #include "nbnxn_ocl_types.h"
64
65 /*! \brief Stringifies the input argument
66  */
67 #define STRINGIFY_PARAM(c) #c
68
69 /*! \brief Stringifies the result of expansion of a macro argument
70  */
71 #define STRINGIFY_MACRO(c) STRINGIFY_PARAM(c)
72
73 /*! \brief Array of the defines needed to generate a specific eel flavour
74  *
75  * The twin-cutoff entries are not normally used, because those setups are
76  * not available to the user. FastGen takes care of generating both
77  * single- and twin-cutoff versions because PME tuning might need both.
78  */
79 static const char * kernel_electrostatic_family_definitions[] =
80 {
81     " -DEL_CUTOFF -DEELNAME=_ElecCut",
82     " -DEL_RF -DEELNAME=_ElecRF",
83     " -DEL_EWALD_TAB -DEELNAME=_ElecEwQSTab",
84     " -DEL_EWALD_TAB -DVDW_CUTOFF_CHECK -DEELNAME=_ElecEwQSTabTwinCut",
85     " -DEL_EWALD_ANA -DEELNAME=_ElecEw",
86     " -DEL_EWALD_ANA -DVDW_CUTOFF_CHECK -DEELNAME=_ElecEwTwinCut"
87 };
88
89 /*! \brief Array of the defines needed to generate a specific vdw flavour
90  */
91 static const char * kernel_VdW_family_definitions[] =
92 {
93     " -DVDWNAME=_VdwLJ",
94     " -DLJ_COMB_GEOM -DVDWNAME=_VdwLJCombGeom",
95     " -DLJ_COMB_LB -DVDWNAME=_VdwLJCombLB",
96     " -DLJ_FORCE_SWITCH -DVDWNAME=_VdwLJFsw",
97     " -DLJ_POT_SWITCH -DVDWNAME=_VdwLJPsw",
98     " -DLJ_EWALD_COMB_GEOM -DVDWNAME=_VdwLJEwCombGeom",
99     " -DLJ_EWALD_COMB_LB -DVDWNAME=_VdwLJEwCombLB"
100 };
101
102 /*! \brief Returns a string with the compiler defines required to avoid all flavour generation
103  *
104  * For example if flavour eelOclRF with evdwOclFSWITCH, the output will be such that the corresponding
105  * kernel flavour is generated:
106  * -DGMX_OCL_FASTGEN          (will replace flavour generator nbnxn_ocl_kernels.clh with nbnxn_ocl_kernels_fastgen.clh)
107  * -DEL_RF                    (The eelOclRF flavour)
108  * -DEELNAME=_ElecRF          (The first part of the generated kernel name )
109  * -DLJ_EWALD_COMB_GEOM       (The evdwOclFSWITCH flavour)
110  * -DVDWNAME=_VdwLJEwCombGeom (The second part of the generated kernel name )
111  *
112  * prune/energy are still generated as originally. It is only the the flavour-level that has changed, so that
113  * only the required flavour for the simulation is compiled.
114  *
115  * If eeltype is single-range Ewald, then we need to add the
116  * twin-cutoff flavour kernels to the JIT, because PME tuning might
117  * need it. This path sets -DGMX_OCL_FASTGEN_ADD_TWINCUT, which
118  * triggers the use of nbnxn_ocl_kernels_fastgen_add_twincut.clh. This
119  * hard-codes the generation of extra kernels that have the same base
120  * flavour, and add the required -DVDW_CUTOFF_CHECK and "TwinCut" to
121  * the kernel name.
122  *
123  * If FastGen is not active, then nothing needs to be returned. The
124  * JIT defaults to compiling all kernel flavours.
125  *
126  * \param[in]  bFastGen    Whether FastGen should be used
127  * \param[in]  eeltype     Electrostatics kernel flavour for FastGen
128  * \param[in]  vdwtype     VDW kernel flavour for FastGen
129  * \return                 String with the defines if FastGen is active
130  *
131  * \throws std::bad_alloc if out of memory
132  */
133 static std::string
134 makeDefinesForKernelTypes(bool bFastGen,
135                           int  eeltype,
136                           int  vdwtype)
137 {
138     std::string defines_for_kernel_types;
139
140     if (bFastGen)
141     {
142         bool bIsEwaldSingleCutoff = (eeltype == eelOclEWALD_TAB ||
143                                      eeltype == eelOclEWALD_ANA);
144
145         if (bIsEwaldSingleCutoff)
146         {
147             defines_for_kernel_types += "-DGMX_OCL_FASTGEN_ADD_TWINCUT";
148         }
149         else
150         {
151             /* This triggers the use of
152                nbnxn_ocl_kernels_fastgen.clh. */
153             defines_for_kernel_types += "-DGMX_OCL_FASTGEN";
154         }
155         defines_for_kernel_types += kernel_electrostatic_family_definitions[eeltype];
156         defines_for_kernel_types += kernel_VdW_family_definitions[vdwtype];
157     }
158
159     return defines_for_kernel_types;
160 }
161
162 /*! \brief Compiles nbnxn kernels for OpenCL GPU given by \p mygpu
163  *
164  * With OpenCL, a call to this function must precede nbnxn_gpu_init().
165  *
166  * Doing bFastGen means only the requested kernels are compiled,
167  * significantly reducing the total compilation time. If false, all
168  * OpenCL kernels are compiled.
169  *
170  * A fatal error results if compilation fails.
171  *
172  * \param[inout] nb  Manages OpenCL non-bonded calculations; compiled kernels returned in dev_info members
173  *
174  * Does not throw
175  */
176 void
177 nbnxn_gpu_compile_kernels(gmx_nbnxn_ocl_t *nb)
178 {
179     gmx_bool                  bFastGen = TRUE;
180     cl_program                program  = nullptr;
181
182     if (getenv("GMX_OCL_NOFASTGEN") != NULL)
183     {
184         bFastGen = FALSE;
185     }
186
187     /* Need to catch std::bad_alloc here and during compilation string
188        handling. */
189     try
190     {
191         std::string extraDefines = makeDefinesForKernelTypes(bFastGen,
192                                                              nb->nbparam->eeltype,
193                                                              nb->nbparam->vdwtype);
194
195         /* Here we pass macros and static const int variables defined in include
196          * files outside the nbnxn_ocl as macros, to avoid including those files
197          * in the JIT compilation that happens at runtime.
198          */
199
200         extraDefines += gmx::formatString(
201                     " -DCENTRAL=%d "
202                     "-DNBNXN_GPU_NCLUSTER_PER_SUPERCLUSTER=%d -DNBNXN_GPU_CLUSTER_SIZE=%d -DNBNXN_GPU_JGROUP_SIZE=%d "
203                     "-DGMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY=%d "
204                     "-DNBNXN_MIN_RSQ=%s %s",
205                     CENTRAL,                                                /* Defined in ishift.h */
206                     c_nbnxnGpuNumClusterPerSupercluster,                    /* Defined in nbnxn_pairlist.h */
207                     c_nbnxnGpuClusterSize,                                  /* Defined in nbnxn_pairlist.h */
208                     c_nbnxnGpuJgroupSize,                                   /* Defined in nbnxn_pairlist.h */
209                     getOclPruneKernelJ4Concurrency(nb->dev_info->vendor_e), /* In nbnxn_ocl_types.h  */
210                     STRINGIFY_MACRO(NBNXN_MIN_RSQ)                          /* Defined in nbnxn_consts.h */
211                                                                             /* NBNXN_MIN_RSQ passed as string to avoid
212                                                                                 floating point representation problems with sprintf */
213                     , (nb->bPrefetchLjParam) ? "-DIATYPE_SHMEM" : ""
214                     );
215
216         try
217         {
218             /* TODO when we have a proper MPI-aware logging module,
219                the log output here should be written there */
220             program = gmx::ocl::compileProgram(stderr,
221                                                "nbnxn_ocl_kernels.cl",
222                                                extraDefines,
223                                                nb->dev_rundata->context,
224                                                nb->dev_info->ocl_gpu_id.ocl_device_id,
225                                                nb->dev_info->vendor_e);
226         }
227         catch (gmx::GromacsException &e)
228         {
229             e.prependContext(gmx::formatString("Failed to compile NBNXN kernels for GPU #%s\n",
230                                                nb->dev_info->device_name));
231             throw;
232         }
233     }
234     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
235
236     nb->dev_rundata->program = program;
237 }