PME spline+spread CUDA kernel and unit tests
[alexxy/gromacs.git] / src / gromacs / ewald / pme-gpu.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017, 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
36 /*! \internal \file
37  * \brief Implements high-level PME GPU functions which do not require GPU framework-specific code.
38  *
39  * \author Aleksei Iupinov <a.yupinov@gmail.com>
40  * \ingroup module_ewald
41  */
42
43 #include "gmxpre.h"
44
45 #include "config.h"
46
47 #include <list>
48
49 #include "gromacs/ewald/pme.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/utility/exceptions.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/stringutil.h"
55
56 #include "pme-gpu-internal.h"
57 #include "pme-grid.h"
58 #include "pme-internal.h"
59 #include "pme-solve.h"
60
61 bool pme_gpu_task_enabled(const gmx_pme_t *pme)
62 {
63     return (pme != nullptr) && (pme->runMode != PmeRunMode::CPU);
64 }
65
66 bool pme_gpu_supports_input(const t_inputrec *ir, std::string *error)
67 {
68     std::list<std::string> errorReasons;
69     if (!EEL_PME(ir->coulombtype))
70     {
71         errorReasons.push_back("systems that do not use PME for electrostatics");
72     }
73     if (ir->pme_order != 4)
74     {
75         errorReasons.push_back("interpolation orders other than 4");
76     }
77     if (ir->efep != efepNO)
78     {
79         errorReasons.push_back("free energy calculations (multiple grids)");
80     }
81     if (EVDW_PME(ir->vdwtype))
82     {
83         errorReasons.push_back("Lennard-Jones PME");
84     }
85 #if GMX_DOUBLE
86     {
87         errorReasons.push_back("double precision");
88     }
89 #endif
90 #if GMX_GPU != GMX_GPU_CUDA
91     {
92         errorReasons.push_back("non-CUDA build of GROMACS");
93     }
94 #endif
95     if (ir->cutoff_scheme == ecutsGROUP)
96     {
97         errorReasons.push_back("group cutoff scheme");
98     }
99     if (EI_TPI(ir->eI))
100     {
101         errorReasons.push_back("test particle insertion");
102     }
103
104     bool inputSupported = errorReasons.empty();
105     if (!inputSupported && error)
106     {
107         std::string regressionTestMarker = "PME GPU does not support";
108         // this prefix is tested for in the regression tests script gmxtest.pl
109         *error = regressionTestMarker + ": " + gmx::joinStrings(errorReasons, "; ") + ".";
110     }
111     return inputSupported;
112 }
113
114 void pme_gpu_reset_timings(const gmx_pme_t *pme)
115 {
116     if (pme_gpu_active(pme))
117     {
118         pme_gpu_reset_timings(pme->gpu);
119     }
120 }
121
122 void pme_gpu_get_timings(const gmx_pme_t *pme, gmx_wallclock_gpu_pme_t *timings)
123 {
124     if (pme_gpu_active(pme))
125     {
126         pme_gpu_get_timings(pme->gpu, timings);
127     }
128 }
129
130 void pme_gpu_get_results(const gmx_pme_t *pme,
131                          gmx_wallcycle_t  wcycle,
132                          matrix           vir_q,
133                          real            *energy_q)
134 {
135     GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
136
137     const bool haveComputedEnergyAndVirial = pme->gpu->settings.stepFlags & GMX_PME_CALC_ENER_VIR;
138     const bool haveComputedForces          = pme->gpu->settings.stepFlags & GMX_PME_CALC_F;
139
140     // The wallcycle hierarchy is messy - we pause ewcFORCE just to exclude the PME GPU sync time
141     wallcycle_stop(wcycle, ewcFORCE);
142
143     wallcycle_start(wcycle, ewcWAIT_GPU_PME_GATHER);
144     pme_gpu_finish_step(pme->gpu, haveComputedForces, haveComputedEnergyAndVirial);
145     wallcycle_stop(wcycle, ewcWAIT_GPU_PME_GATHER);
146
147     wallcycle_start_nocount(wcycle, ewcFORCE);
148
149     if (haveComputedEnergyAndVirial)
150     {
151         if (pme->doCoulomb)
152         {
153             if (pme_gpu_performs_solve(pme->gpu))
154             {
155                 pme_gpu_get_energy_virial(pme->gpu, energy_q, vir_q);
156             }
157             else
158             {
159                 get_pme_ener_vir_q(pme->solve_work, pme->nthread, energy_q, vir_q);
160             }
161         }
162         else
163         {
164             *energy_q = 0;
165         }
166     }
167     /* No additional haveComputedForces code since forces are copied to the output host buffer with no transformation. */
168 }