Set up workload data structures
[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,2018,2019, 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/ewald_utils.h"
50 #include "gromacs/ewald/pme.h"
51 #include "gromacs/fft/parallel_3dfft.h"
52 #include "gromacs/math/invertmatrix.h"
53 #include "gromacs/mdlib/gmx_omp_nthreads.h"
54 #include "gromacs/mdtypes/enerdata.h"
55 #include "gromacs/mdtypes/forceoutput.h"
56 #include "gromacs/mdtypes/inputrec.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/stringutil.h"
61
62 #include "pme_gpu_internal.h"
63 #include "pme_grid.h"
64 #include "pme_internal.h"
65 #include "pme_solve.h"
66
67 void pme_gpu_reset_timings(const gmx_pme_t *pme)
68 {
69     if (pme_gpu_active(pme))
70     {
71         pme_gpu_reset_timings(pme->gpu);
72     }
73 }
74
75 void pme_gpu_get_timings(const gmx_pme_t *pme, gmx_wallclock_gpu_pme_t *timings)
76 {
77     if (pme_gpu_active(pme))
78     {
79         pme_gpu_get_timings(pme->gpu, timings);
80     }
81 }
82
83 int pme_gpu_get_padding_size(const gmx_pme_t *pme)
84 {
85
86     if (!pme || !pme_gpu_active(pme))
87     {
88         return 0;
89     }
90     else
91     {
92         return pme_gpu_get_atom_data_alignment(pme->gpu);
93     }
94 }
95
96 /*! \brief
97  * A convenience wrapper for launching either the GPU or CPU FFT.
98  *
99  * \param[in] pme            The PME structure.
100  * \param[in] gridIndex      The grid index - should currently always be 0.
101  * \param[in] dir            The FFT direction enum.
102  * \param[in] wcycle         The wallclock counter.
103  */
104 void inline parallel_3dfft_execute_gpu_wrapper(gmx_pme_t              *pme,
105                                                const int               gridIndex,
106                                                enum gmx_fft_direction  dir,
107                                                gmx_wallcycle_t         wcycle)
108 {
109     GMX_ASSERT(gridIndex == 0, "Only single grid supported");
110     if (pme_gpu_performs_FFT(pme->gpu))
111     {
112         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
113         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
114         pme_gpu_3dfft(pme->gpu, dir, gridIndex);
115         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
116         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
117     }
118     else
119     {
120         wallcycle_start(wcycle, ewcPME_FFT_MIXED_MODE);
121 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
122         for (int thread = 0; thread < pme->nthread; thread++)
123         {
124             gmx_parallel_3dfft_execute(pme->pfft_setup[gridIndex], dir, thread, wcycle);
125         }
126         wallcycle_stop(wcycle, ewcPME_FFT_MIXED_MODE);
127     }
128 }
129
130 /* The PME computation code split into a few separate functions. */
131
132 void pme_gpu_prepare_computation(gmx_pme_t            *pme,
133                                  bool                  needToUpdateBox,
134                                  const matrix          box,
135                                  gmx_wallcycle        *wcycle,
136                                  int                   flags,
137                                  bool                  useGpuForceReduction)
138 {
139     GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
140     GMX_ASSERT(pme->nnodes > 0, "");
141     GMX_ASSERT(pme->nnodes == 1 || pme->ndecompdim > 0, "");
142
143     PmeGpu *pmeGpu = pme->gpu;
144     pmeGpu->settings.currentFlags         = flags;
145     // TODO these flags are only here to honor the CPU PME code, and probably should be removed
146     pmeGpu->settings.useGpuForceReduction = useGpuForceReduction;
147
148     bool shouldUpdateBox = false;
149     for (int i = 0; i < DIM; ++i)
150     {
151         for (int j = 0; j <= i; ++j)
152         {
153             shouldUpdateBox                  |= (pmeGpu->common->previousBox[i][j] != box[i][j]);
154             pmeGpu->common->previousBox[i][j] = box[i][j];
155         }
156     }
157
158     if (needToUpdateBox || shouldUpdateBox) // || is to make the first computation always update
159     {
160         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
161         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
162         pme_gpu_update_input_box(pmeGpu, box);
163         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
164         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
165
166         if (!pme_gpu_performs_solve(pmeGpu))
167         {
168             // TODO remove code duplication and add test coverage
169             matrix scaledBox;
170             pmeGpu->common->boxScaler->scaleBox(box, scaledBox);
171             gmx::invertBoxMatrix(scaledBox, pme->recipbox);
172             pme->boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
173         }
174     }
175 }
176
177 void pme_gpu_copy_coordinates_to_gpu(gmx_pme_t            *pme,
178                                      const rvec           *coordinatesHost,
179                                      gmx_wallcycle        *wcycle)
180 {
181     GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
182
183     PmeGpu *pmeGpu = pme->gpu;
184
185     // The only spot of PME GPU where LAUNCH_GPU counter increases call-count
186     wallcycle_start(wcycle, ewcLAUNCH_GPU);
187     // The only spot of PME GPU where ewcsLAUNCH_GPU_PME subcounter increases call-count
188     wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_PME);
189     pme_gpu_copy_input_coordinates(pmeGpu, coordinatesHost);
190     wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
191     wallcycle_stop(wcycle, ewcLAUNCH_GPU);
192 }
193
194 void pme_gpu_launch_spread(gmx_pme_t            *pme,
195                            gmx_wallcycle        *wcycle)
196 {
197     GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
198
199     PmeGpu            *pmeGpu = pme->gpu;
200
201     const unsigned int gridIndex  = 0;
202     real              *fftgrid    = pme->fftgrid[gridIndex];
203     if (pmeGpu->settings.currentFlags & GMX_PME_SPREAD)
204     {
205         /* Spread the coefficients on a grid */
206         const bool computeSplines = true;
207         const bool spreadCharges  = true;
208         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
209         wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
210         pme_gpu_spread(pmeGpu, gridIndex, fftgrid, computeSplines, spreadCharges);
211         wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
212         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
213     }
214 }
215
216 void pme_gpu_launch_complex_transforms(gmx_pme_t      *pme,
217                                        gmx_wallcycle  *wcycle)
218 {
219     PmeGpu            *pmeGpu                 = pme->gpu;
220     const bool         computeEnergyAndVirial = (pmeGpu->settings.currentFlags & GMX_PME_CALC_ENER_VIR) != 0;
221     const bool         performBackFFT         = (pmeGpu->settings.currentFlags & (GMX_PME_CALC_F | GMX_PME_CALC_POT)) != 0;
222     const unsigned int gridIndex              = 0;
223     t_complex         *cfftgrid               = pme->cfftgrid[gridIndex];
224
225     if (pmeGpu->settings.currentFlags & GMX_PME_SPREAD)
226     {
227         if (!pme_gpu_performs_FFT(pmeGpu))
228         {
229             wallcycle_start(wcycle, ewcWAIT_GPU_PME_SPREAD);
230             pme_gpu_sync_spread_grid(pme->gpu);
231             wallcycle_stop(wcycle, ewcWAIT_GPU_PME_SPREAD);
232         }
233     }
234
235     try
236     {
237         if (pmeGpu->settings.currentFlags & GMX_PME_SOLVE)
238         {
239             /* do R2C 3D-FFT */
240             parallel_3dfft_execute_gpu_wrapper(pme, gridIndex, GMX_FFT_REAL_TO_COMPLEX, wcycle);
241
242             /* solve in k-space for our local cells */
243             if (pme_gpu_performs_solve(pmeGpu))
244             {
245                 const auto gridOrdering = pme_gpu_uses_dd(pmeGpu) ? GridOrdering::YZX : GridOrdering::XYZ;
246                 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
247                 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
248                 pme_gpu_solve(pmeGpu, cfftgrid, gridOrdering, computeEnergyAndVirial);
249                 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
250                 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
251             }
252             else
253             {
254                 wallcycle_start(wcycle, ewcPME_SOLVE_MIXED_MODE);
255 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
256                 for (int thread = 0; thread < pme->nthread; thread++)
257                 {
258                     solve_pme_yzx(pme, cfftgrid, pme->boxVolume,
259                                   computeEnergyAndVirial, pme->nthread, thread);
260                 }
261                 wallcycle_stop(wcycle, ewcPME_SOLVE_MIXED_MODE);
262             }
263         }
264
265         if (performBackFFT)
266         {
267             parallel_3dfft_execute_gpu_wrapper(pme, gridIndex, GMX_FFT_COMPLEX_TO_REAL, wcycle);
268         }
269     } GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
270 }
271
272 void pme_gpu_launch_gather(const gmx_pme_t                 *pme,
273                            gmx_wallcycle gmx_unused        *wcycle,
274                            PmeForceOutputHandling           forceTreatment)
275 {
276     GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
277
278     if (!pme_gpu_performs_gather(pme->gpu))
279     {
280         return;
281     }
282
283     wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
284     wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
285     const unsigned int gridIndex  = 0;
286     real              *fftgrid    = pme->fftgrid[gridIndex];
287     pme_gpu_gather(pme->gpu, forceTreatment, reinterpret_cast<float *>(fftgrid));
288     wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
289     wallcycle_stop(wcycle, ewcLAUNCH_GPU);
290 }
291
292 //! Accumulate the \c forcesToAdd to \c f, using the available threads.
293 static void sum_forces(gmx::ArrayRef<gmx::RVec>       f,
294                        gmx::ArrayRef<const gmx::RVec> forceToAdd)
295 {
296     const int      end = forceToAdd.size();
297
298     int gmx_unused nt = gmx_omp_nthreads_get(emntPME);
299 #pragma omp parallel for num_threads(nt) schedule(static)
300     for (int i = 0; i < end; i++)
301     {
302         f[i] += forceToAdd[i];
303     }
304 }
305
306 //! Reduce quantities from \c output to \c forceWithVirial and \c enerd.
307 static void pme_gpu_reduce_outputs(const int             flags,
308                                    const PmeOutput      &output,
309                                    gmx_wallcycle        *wcycle,
310                                    gmx::ForceWithVirial *forceWithVirial,
311                                    gmx_enerdata_t       *enerd)
312 {
313     wallcycle_start(wcycle, ewcPME_GPU_F_REDUCTION);
314     GMX_ASSERT(forceWithVirial, "Invalid force pointer");
315
316     const bool haveComputedEnergyAndVirial = (flags & GMX_PME_CALC_ENER_VIR) != 0;
317     if (haveComputedEnergyAndVirial)
318     {
319         GMX_ASSERT(enerd, "Invalid energy output manager");
320         forceWithVirial->addVirialContribution(output.coulombVirial_);
321         enerd->term[F_COUL_RECIP] += output.coulombEnergy_;
322     }
323     if (output.haveForceOutput_)
324     {
325         sum_forces(forceWithVirial->force_, output.forces_);
326     }
327     wallcycle_stop(wcycle, ewcPME_GPU_F_REDUCTION);
328 }
329
330 bool pme_gpu_try_finish_task(gmx_pme_t            *pme,
331                              const int             flags,
332                              gmx_wallcycle        *wcycle,
333                              gmx::ForceWithVirial *forceWithVirial,
334                              gmx_enerdata_t       *enerd,
335                              GpuTaskCompletion     completionKind)
336 {
337     GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
338     GMX_ASSERT(!pme->gpu->settings.useGpuForceReduction, "GPU force reduction should not be active on the pme_gpu_try_finish_task() path");
339
340     // First, if possible, check whether all tasks on the stream have
341     // completed, and return fast if not. Accumulate to wcycle the
342     // time needed for that checking, but do not yet record that the
343     // gather has occured.
344     bool           needToSynchronize      = true;
345     constexpr bool c_streamQuerySupported = (GMX_GPU == GMX_GPU_CUDA);
346     // TODO: implement c_streamQuerySupported with an additional GpuEventSynchronizer per stream (#2521)
347     if ((completionKind == GpuTaskCompletion::Check) && c_streamQuerySupported)
348     {
349         wallcycle_start_nocount(wcycle, ewcWAIT_GPU_PME_GATHER);
350         // Query the PME stream for completion of all tasks enqueued and
351         // if we're not done, stop the timer before early return.
352         const bool pmeGpuDone = pme_gpu_stream_query(pme->gpu);
353         wallcycle_stop(wcycle, ewcWAIT_GPU_PME_GATHER);
354
355         if (!pmeGpuDone)
356         {
357             return false;
358         }
359         needToSynchronize = false;
360     }
361
362     wallcycle_start(wcycle, ewcWAIT_GPU_PME_GATHER);
363     // If the above check passed, then there is no need to make an
364     // explicit synchronization call.
365     if (needToSynchronize)
366     {
367         // Synchronize the whole PME stream at once, including D2H result transfers.
368         pme_gpu_synchronize(pme->gpu);
369     }
370     pme_gpu_update_timings(pme->gpu);
371     PmeOutput output = pme_gpu_getOutput(*pme, flags);
372     wallcycle_stop(wcycle, ewcWAIT_GPU_PME_GATHER);
373
374     GMX_ASSERT(pme->gpu->settings.useGpuForceReduction == !output.haveForceOutput_, "When forces are reduced on the CPU, there needs to be force output");
375     pme_gpu_reduce_outputs(flags, output, wcycle, forceWithVirial, enerd);
376
377     return true;
378 }
379
380 // This is used by PME-only ranks
381 PmeOutput pme_gpu_wait_finish_task(gmx_pme_t     *pme,
382                                    const int      flags,
383                                    gmx_wallcycle *wcycle)
384 {
385     GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
386
387     wallcycle_start(wcycle, ewcWAIT_GPU_PME_GATHER);
388
389     // Synchronize the whole PME stream at once, including D2H result transfers.
390     //
391     // TODO: make this sync conditional to wait only for virial/energies
392     pme_gpu_synchronize(pme->gpu);
393
394     pme_gpu_update_timings(pme->gpu);
395     PmeOutput output = pme_gpu_getOutput(*pme, flags);
396     wallcycle_stop(wcycle, ewcWAIT_GPU_PME_GATHER);
397     return output;
398 }
399
400 // This is used when not using the alternate-waiting reduction
401 void pme_gpu_wait_and_reduce(gmx_pme_t            *pme,
402                              const int             flags,
403                              gmx_wallcycle        *wcycle,
404                              gmx::ForceWithVirial *forceWithVirial,
405                              gmx_enerdata_t       *enerd)
406 {
407     PmeOutput output = pme_gpu_wait_finish_task(pme, flags, wcycle);
408     GMX_ASSERT(pme->gpu->settings.useGpuForceReduction == !output.haveForceOutput_, "When forces are reduced on the CPU, there needs to be force output");
409     pme_gpu_reduce_outputs(flags, output, wcycle, forceWithVirial, enerd);
410 }
411
412 void pme_gpu_reinit_computation(const gmx_pme_t *pme,
413                                 gmx_wallcycle   *wcycle)
414 {
415     GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
416
417     wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
418     wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
419
420     pme_gpu_clear_grids(pme->gpu);
421     pme_gpu_clear_energy_virial(pme->gpu);
422
423     wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
424     wallcycle_stop(wcycle, ewcLAUNCH_GPU);
425 }
426
427 DeviceBuffer<float> pme_gpu_get_device_x(const gmx_pme_t *pme)
428 {
429     GMX_ASSERT((pme && pme_gpu_active(pme)), "PME GPU coordinates buffer was requested from uninitialized PME module");
430     return pme_gpu_get_kernelparam_coordinates(pme->gpu);
431 }
432
433 void *pme_gpu_get_device_f(const gmx_pme_t *pme)
434 {
435     if (!pme || !pme_gpu_active(pme))
436     {
437         return nullptr;
438     }
439     return pme_gpu_get_kernelparam_forces(pme->gpu);
440 }
441
442 void *pme_gpu_get_device_stream(const gmx_pme_t *pme)
443 {
444     if (!pme || !pme_gpu_active(pme))
445     {
446         return nullptr;
447     }
448     return pme_gpu_get_stream(pme->gpu);
449 }
450
451 GpuEventSynchronizer * pme_gpu_get_f_ready_synchronizer(const gmx_pme_t *pme)
452 {
453     if (!pme || !pme_gpu_active(pme))
454     {
455         return nullptr;
456     }
457
458     return pme_gpu_get_forces_ready_synchronizer(pme->gpu);
459 }