2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * \brief Implements high-level PME GPU functions which do not require GPU framework-specific code.
39 * \author Aleksei Iupinov <a.yupinov@gmail.com>
40 * \ingroup module_ewald
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"
62 #include "pme_gpu_internal.h"
64 #include "pme_internal.h"
65 #include "pme_solve.h"
67 void pme_gpu_reset_timings(const gmx_pme_t *pme)
69 if (pme_gpu_active(pme))
71 pme_gpu_reset_timings(pme->gpu);
75 void pme_gpu_get_timings(const gmx_pme_t *pme, gmx_wallclock_gpu_pme_t *timings)
77 if (pme_gpu_active(pme))
79 pme_gpu_get_timings(pme->gpu, timings);
83 int pme_gpu_get_padding_size(const gmx_pme_t *pme)
86 if (!pme || !pme_gpu_active(pme))
92 return pme_gpu_get_atom_data_alignment(pme->gpu);
97 * A convenience wrapper for launching either the GPU or CPU FFT.
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.
104 void inline parallel_3dfft_execute_gpu_wrapper(gmx_pme_t *pme,
106 enum gmx_fft_direction dir,
107 gmx_wallcycle_t wcycle)
109 GMX_ASSERT(gridIndex == 0, "Only single grid supported");
110 if (pme_gpu_performs_FFT(pme->gpu))
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);
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++)
124 gmx_parallel_3dfft_execute(pme->pfft_setup[gridIndex], dir, thread, wcycle);
126 wallcycle_stop(wcycle, ewcPME_FFT_MIXED_MODE);
130 /* The PME computation code split into a few separate functions. */
132 void pme_gpu_prepare_computation(gmx_pme_t *pme,
133 bool needToUpdateBox,
135 gmx_wallcycle *wcycle,
138 GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
139 GMX_ASSERT(pme->nnodes > 0, "");
140 GMX_ASSERT(pme->nnodes == 1 || pme->ndecompdim > 0, "");
142 PmeGpu *pmeGpu = pme->gpu;
143 pmeGpu->settings.currentFlags = flags;
144 // TODO these flags are only here to honor the CPU PME code, and probably should be removed
146 bool shouldUpdateBox = false;
147 for (int i = 0; i < DIM; ++i)
149 for (int j = 0; j <= i; ++j)
151 shouldUpdateBox |= (pmeGpu->common->previousBox[i][j] != box[i][j]);
152 pmeGpu->common->previousBox[i][j] = box[i][j];
156 if (needToUpdateBox || shouldUpdateBox) // || is to make the first computation always update
158 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
159 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
160 pme_gpu_update_input_box(pmeGpu, box);
161 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
162 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
164 if (!pme_gpu_performs_solve(pmeGpu))
166 // TODO remove code duplication and add test coverage
168 pmeGpu->common->boxScaler->scaleBox(box, scaledBox);
169 gmx::invertBoxMatrix(scaledBox, pme->recipbox);
170 pme->boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
175 void pme_gpu_copy_coordinates_to_gpu(gmx_pme_t *pme,
176 const rvec *coordinatesHost,
177 gmx_wallcycle *wcycle)
179 GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
181 PmeGpu *pmeGpu = pme->gpu;
183 // The only spot of PME GPU where LAUNCH_GPU counter increases call-count
184 wallcycle_start(wcycle, ewcLAUNCH_GPU);
185 // The only spot of PME GPU where ewcsLAUNCH_GPU_PME subcounter increases call-count
186 wallcycle_sub_start(wcycle, ewcsLAUNCH_GPU_PME);
187 pme_gpu_copy_input_coordinates(pmeGpu, coordinatesHost);
188 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
189 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
192 void pme_gpu_launch_spread(gmx_pme_t *pme,
193 gmx_wallcycle *wcycle)
195 GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
197 PmeGpu *pmeGpu = pme->gpu;
199 const unsigned int gridIndex = 0;
200 real *fftgrid = pme->fftgrid[gridIndex];
201 if (pmeGpu->settings.currentFlags & GMX_PME_SPREAD)
203 /* Spread the coefficients on a grid */
204 const bool computeSplines = true;
205 const bool spreadCharges = true;
206 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
207 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
208 pme_gpu_spread(pmeGpu, gridIndex, fftgrid, computeSplines, spreadCharges);
209 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
210 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
214 void pme_gpu_launch_complex_transforms(gmx_pme_t *pme,
215 gmx_wallcycle *wcycle)
217 PmeGpu *pmeGpu = pme->gpu;
218 const bool computeEnergyAndVirial = (pmeGpu->settings.currentFlags & GMX_PME_CALC_ENER_VIR) != 0;
219 const bool performBackFFT = (pmeGpu->settings.currentFlags & (GMX_PME_CALC_F | GMX_PME_CALC_POT)) != 0;
220 const unsigned int gridIndex = 0;
221 t_complex *cfftgrid = pme->cfftgrid[gridIndex];
223 if (pmeGpu->settings.currentFlags & GMX_PME_SPREAD)
225 if (!pme_gpu_performs_FFT(pmeGpu))
227 wallcycle_start(wcycle, ewcWAIT_GPU_PME_SPREAD);
228 pme_gpu_sync_spread_grid(pme->gpu);
229 wallcycle_stop(wcycle, ewcWAIT_GPU_PME_SPREAD);
235 if (pmeGpu->settings.currentFlags & GMX_PME_SOLVE)
238 parallel_3dfft_execute_gpu_wrapper(pme, gridIndex, GMX_FFT_REAL_TO_COMPLEX, wcycle);
240 /* solve in k-space for our local cells */
241 if (pme_gpu_performs_solve(pmeGpu))
243 const auto gridOrdering = pme_gpu_uses_dd(pmeGpu) ? GridOrdering::YZX : GridOrdering::XYZ;
244 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
245 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
246 pme_gpu_solve(pmeGpu, cfftgrid, gridOrdering, computeEnergyAndVirial);
247 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
248 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
252 wallcycle_start(wcycle, ewcPME_SOLVE_MIXED_MODE);
253 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
254 for (int thread = 0; thread < pme->nthread; thread++)
256 solve_pme_yzx(pme, cfftgrid, pme->boxVolume,
257 computeEnergyAndVirial, pme->nthread, thread);
259 wallcycle_stop(wcycle, ewcPME_SOLVE_MIXED_MODE);
265 parallel_3dfft_execute_gpu_wrapper(pme, gridIndex, GMX_FFT_COMPLEX_TO_REAL, wcycle);
267 } GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
270 void pme_gpu_launch_gather(const gmx_pme_t *pme,
271 gmx_wallcycle gmx_unused *wcycle,
272 PmeForceOutputHandling forceTreatment,
273 bool useGpuFPmeReduction)
275 GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
277 if (!pme_gpu_performs_gather(pme->gpu))
282 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
283 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
284 const unsigned int gridIndex = 0;
285 real *fftgrid = pme->fftgrid[gridIndex];
286 pme_gpu_gather(pme->gpu, forceTreatment, reinterpret_cast<float *>(fftgrid), useGpuFPmeReduction);
287 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
288 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
291 //! Accumulate the \c forcesToAdd to \c f, using the available threads.
292 static void sum_forces(gmx::ArrayRef<gmx::RVec> f,
293 gmx::ArrayRef<const gmx::RVec> forceToAdd)
295 const int end = forceToAdd.size();
297 int gmx_unused nt = gmx_omp_nthreads_get(emntPME);
298 #pragma omp parallel for num_threads(nt) schedule(static)
299 for (int i = 0; i < end; i++)
301 f[i] += forceToAdd[i];
305 //! Reduce quantities from \c output to \c forceWithVirial and \c enerd.
306 static void pme_gpu_reduce_outputs(const int flags,
307 const PmeOutput &output,
308 gmx_wallcycle *wcycle,
309 gmx::ForceWithVirial *forceWithVirial,
310 gmx_enerdata_t *enerd,
311 bool useGpuFPmeReduction)
313 wallcycle_start(wcycle, ewcPME_GPU_F_REDUCTION);
314 GMX_ASSERT(forceWithVirial, "Invalid force pointer");
315 const bool haveComputedEnergyAndVirial = (flags & GMX_PME_CALC_ENER_VIR) != 0;
316 if (haveComputedEnergyAndVirial)
318 GMX_ASSERT(enerd, "Invalid energy output manager");
319 forceWithVirial->addVirialContribution(output.coulombVirial_);
320 enerd->term[F_COUL_RECIP] += output.coulombEnergy_;
322 if (!useGpuFPmeReduction)
324 sum_forces(forceWithVirial->force_, output.forces_);
326 wallcycle_stop(wcycle, ewcPME_GPU_F_REDUCTION);
329 bool pme_gpu_try_finish_task(gmx_pme_t *pme,
331 gmx_wallcycle *wcycle,
332 gmx::ForceWithVirial *forceWithVirial,
333 gmx_enerdata_t *enerd,
334 GpuTaskCompletion completionKind)
336 GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
338 // First, if possible, check whether all tasks on the stream have
339 // completed, and return fast if not. Accumulate to wcycle the
340 // time needed for that checking, but do not yet record that the
341 // gather has occured.
342 bool needToSynchronize = true;
343 constexpr bool c_streamQuerySupported = (GMX_GPU == GMX_GPU_CUDA);
344 // TODO: implement c_streamQuerySupported with an additional GpuEventSynchronizer per stream (#2521)
345 if ((completionKind == GpuTaskCompletion::Check) && c_streamQuerySupported)
347 wallcycle_start_nocount(wcycle, ewcWAIT_GPU_PME_GATHER);
348 // Query the PME stream for completion of all tasks enqueued and
349 // if we're not done, stop the timer before early return.
350 const bool pmeGpuDone = pme_gpu_stream_query(pme->gpu);
351 wallcycle_stop(wcycle, ewcWAIT_GPU_PME_GATHER);
357 needToSynchronize = false;
360 wallcycle_start(wcycle, ewcWAIT_GPU_PME_GATHER);
361 // If the above check passed, then there is no need to make an
362 // explicit synchronization call.
363 if (needToSynchronize)
365 // Synchronize the whole PME stream at once, including D2H result transfers.
366 pme_gpu_synchronize(pme->gpu);
368 pme_gpu_update_timings(pme->gpu);
369 PmeOutput output = pme_gpu_getOutput(*pme, flags);
370 wallcycle_stop(wcycle, ewcWAIT_GPU_PME_GATHER);
372 pme_gpu_reduce_outputs(flags, output, wcycle, forceWithVirial, enerd, false);
377 // This is used by PME-only ranks
378 PmeOutput pme_gpu_wait_finish_task(gmx_pme_t *pme,
380 gmx_wallcycle *wcycle)
382 GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
384 wallcycle_start(wcycle, ewcWAIT_GPU_PME_GATHER);
386 // Synchronize the whole PME stream at once, including D2H result transfers.
387 // TODO: make this sync conditional with useGpuFPmeReduction to wait only for virial/energies
388 pme_gpu_synchronize(pme->gpu);
390 pme_gpu_update_timings(pme->gpu);
391 PmeOutput output = pme_gpu_getOutput(*pme, flags);
392 wallcycle_stop(wcycle, ewcWAIT_GPU_PME_GATHER);
396 // This is used when not using the alternate-waiting reduction
397 void pme_gpu_wait_and_reduce(gmx_pme_t *pme,
399 gmx_wallcycle *wcycle,
400 gmx::ForceWithVirial *forceWithVirial,
401 gmx_enerdata_t *enerd,
402 bool useGpuFPmeReduction)
404 PmeOutput output = pme_gpu_wait_finish_task(pme, flags, wcycle);
405 pme_gpu_reduce_outputs(flags, output, wcycle, forceWithVirial, enerd, useGpuFPmeReduction);
408 void pme_gpu_reinit_computation(const gmx_pme_t *pme,
409 gmx_wallcycle *wcycle)
411 GMX_ASSERT(pme_gpu_active(pme), "This should be a GPU run of PME but it is not enabled.");
413 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU);
414 wallcycle_sub_start_nocount(wcycle, ewcsLAUNCH_GPU_PME);
416 pme_gpu_clear_grids(pme->gpu);
417 pme_gpu_clear_energy_virial(pme->gpu);
419 wallcycle_sub_stop(wcycle, ewcsLAUNCH_GPU_PME);
420 wallcycle_stop(wcycle, ewcLAUNCH_GPU);
423 DeviceBuffer<float> pme_gpu_get_device_x(const gmx_pme_t *pme)
425 GMX_ASSERT((pme && pme_gpu_active(pme)), "PME GPU coordinates buffer was requested from uninitialized PME module");
426 return pme_gpu_get_kernelparam_coordinates(pme->gpu);
429 void *pme_gpu_get_device_f(const gmx_pme_t *pme)
431 if (!pme || !pme_gpu_active(pme))
435 return pme_gpu_get_kernelparam_forces(pme->gpu);
438 void *pme_gpu_get_device_stream(const gmx_pme_t *pme)
440 if (!pme || !pme_gpu_active(pme))
444 return pme_gpu_get_stream(pme->gpu);
447 GpuEventSynchronizer * pme_gpu_get_f_ready_synchronizer(const gmx_pme_t *pme)
449 if (!pme || !pme_gpu_active(pme))
454 return pme_gpu_get_forces_ready_synchronizer(pme->gpu);