Gromacs  2018.8
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions
#include "gmxpre.h"
#include "config.h"
#include <list>
#include "gromacs/ewald/ewald-utils.h"
#include "gromacs/ewald/pme.h"
#include "gromacs/fft/parallel_3dfft.h"
#include "gromacs/math/invertmatrix.h"
#include "gromacs/mdtypes/inputrec.h"
#include "gromacs/utility/exceptions.h"
#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/gmxassert.h"
#include "gromacs/utility/stringutil.h"
#include "pme-gpu-internal.h"
#include "pme-grid.h"
#include "pme-internal.h"
#include "pme-solve.h"
+ Include dependency graph for pme-gpu.cpp:

Description

Implements high-level PME GPU functions which do not require GPU framework-specific code.

Author
Aleksei Iupinov a.yup.nosp@m.inov.nosp@m.@gmai.nosp@m.l.co.nosp@m.m

Functions

PmeRunMode pme_run_mode (const gmx_pme_t *pme)
 Returns the active PME codepath (CPU, GPU, mixed). More...
 
bool pme_gpu_supports_input (const t_inputrec *ir, std::string *error)
 Checks whether the input system allows to run PME on GPU. TODO: this mostly duplicates an internal PME assert function pme_gpu_check_restrictions(), except that works with a formed gmx_pme_t structure. Should that one go away/work with inputrec? More...
 
void pme_gpu_reset_timings (const gmx_pme_t *pme)
 Resets the PME GPU timings. To be called at the reset step. More...
 
void pme_gpu_get_timings (const gmx_pme_t *pme, gmx_wallclock_gpu_pme_t *timings)
 Copies the PME GPU timings to the gmx_wallclock_gpu_pme_t structure (for log output). To be called at the run end. More...
 
void parallel_3dfft_execute_gpu_wrapper (gmx_pme_t *pme, const int gridIndex, enum gmx_fft_direction dir, gmx_wallcycle_t wcycle)
 A convenience wrapper for launching either the GPU or CPU FFT. More...
 
void pme_gpu_prepare_computation (gmx_pme_t *pme, bool needToUpdateBox, const matrix box, gmx_wallcycle_t wcycle, int flags)
 Prepares PME on GPU computation (updating the box if needed) More...
 
void pme_gpu_launch_spread (gmx_pme_t *pme, const rvec *x, gmx_wallcycle_t wcycle)
 Launches first stage of PME on GPU - H2D input transfers, spreading kernel, and D2H grid transfer if needed. More...
 
void pme_gpu_launch_complex_transforms (gmx_pme_t *pme, gmx_wallcycle_t wcycle)
 Launches middle stages of PME (FFT R2C, solving, FFT C2R) either on GPU or on CPU, depending on the run mode. More...
 
void pme_gpu_launch_gather (const gmx_pme_t *pme, gmx_wallcycle_t wcycle, PmeForceOutputHandling forceTreatment)
 Launches last stage of PME on GPU - force gathering and D2H force transfer. More...
 
static void pme_gpu_get_staged_results (const gmx_pme_t *pme, gmx::ArrayRef< const gmx::RVec > *forces, matrix virial, real *energy)
 Reduce staged virial and energy outputs. More...
 
bool pme_gpu_try_finish_task (const gmx_pme_t *pme, gmx_wallcycle_t wcycle, gmx::ArrayRef< const gmx::RVec > *forces, matrix virial, real *energy, GpuTaskCompletion completionKind)
 Attempts to complete PME GPU tasks. More...
 
void pme_gpu_wait_finish_task (const gmx_pme_t *pme, gmx_wallcycle_t wcycle, gmx::ArrayRef< const gmx::RVec > *forces, matrix virial, real *energy)
 Blocks until PME GPU tasks are completed, and gets the output forces and virial/energy (if they were to be computed). More...
 

Function Documentation

void parallel_3dfft_execute_gpu_wrapper ( gmx_pme_t *  pme,
const int  gridIndex,
enum gmx_fft_direction  dir,
gmx_wallcycle_t  wcycle 
)
inline

A convenience wrapper for launching either the GPU or CPU FFT.

Parameters
[in]pmeThe PME structure.
[in]gridIndexThe grid index - should currently always be 0.
[in]dirThe FFT direction enum.
[in]wcycleThe wallclock counter.
static void pme_gpu_get_staged_results ( const gmx_pme_t *  pme,
gmx::ArrayRef< const gmx::RVec > *  forces,
matrix  virial,
real energy 
)
static

Reduce staged virial and energy outputs.

Parameters
[in]pmeThe PME data structure.
[out]forcesOutput forces pointer, the internal ArrayRef pointers gets assigned to it.
[out]virialThe output virial matrix.
[out]energyThe output energy.
void pme_gpu_get_timings ( const gmx_pme_t *  pme,
gmx_wallclock_gpu_pme_t timings 
)

Copies the PME GPU timings to the gmx_wallclock_gpu_pme_t structure (for log output). To be called at the run end.

Parameters
[in]pmeThe PME structure.
[in]timingsThe gmx_wallclock_gpu_pme_t structure.
void pme_gpu_launch_complex_transforms ( gmx_pme_t *  pme,
gmx_wallcycle_t  wcycle 
)

Launches middle stages of PME (FFT R2C, solving, FFT C2R) either on GPU or on CPU, depending on the run mode.

Parameters
[in]pmeThe PME data structure.
[in]wcycleThe wallclock counter.
void pme_gpu_launch_gather ( const gmx_pme_t *  pme,
gmx_wallcycle_t  wcycle,
PmeForceOutputHandling  forceTreatment 
)

Launches last stage of PME on GPU - force gathering and D2H force transfer.

Parameters
[in]pmeThe PME data structure.
[in]wcycleThe wallclock counter.
[in]forceTreatmentTells how data should be treated. The gathering kernel either stores the output reciprocal forces into the host array, or copies its contents to the GPU first and accumulates. The reduction is non-atomic.
void pme_gpu_launch_spread ( gmx_pme_t *  pme,
const rvec *  x,
gmx_wallcycle_t  wcycle 
)

Launches first stage of PME on GPU - H2D input transfers, spreading kernel, and D2H grid transfer if needed.

Parameters
[in]pmeThe PME data structure.
[in]xThe array of local atoms' coordinates.
[in]wcycleThe wallclock counter.
void pme_gpu_prepare_computation ( gmx_pme_t *  pme,
bool  needToUpdateBox,
const matrix  box,
gmx_wallcycle_t  wcycle,
int  flags 
)

Prepares PME on GPU computation (updating the box if needed)

Parameters
[in]pmeThe PME data structure.
[in]needToUpdateBoxTells if the stored unit cell parameters should be updated from box.
[in]boxThe unit cell box.
[in]wcycleThe wallclock counter.
[in]flagsThe combination of flags to affect this PME computation. The flags are the GMX_PME_ flags from pme.h.
void pme_gpu_reset_timings ( const gmx_pme_t *  pme)

Resets the PME GPU timings. To be called at the reset step.

Parameters
[in]pmeThe PME structure.
bool pme_gpu_supports_input ( const t_inputrec *  ir,
std::string *  error 
)

Checks whether the input system allows to run PME on GPU. TODO: this mostly duplicates an internal PME assert function pme_gpu_check_restrictions(), except that works with a formed gmx_pme_t structure. Should that one go away/work with inputrec?

Parameters
[in]irInput system.
[out]errorThe error message if the input is not supported on GPU.
Returns
true if PME can run on GPU with this input, false otherwise.
bool pme_gpu_try_finish_task ( const gmx_pme_t *  pme,
gmx_wallcycle_t  wcycle,
gmx::ArrayRef< const gmx::RVec > *  forces,
matrix  virial,
real energy,
GpuTaskCompletion  completionKind 
)

Attempts to complete PME GPU tasks.

The completionKind argument controls whether the function blocks until all PME GPU tasks enqueued completed (as pme_gpu_wait_finish_task() does) or only checks and returns immediately if they did not. When blocking or the tasks have completed it also gets the output forces by assigning the ArrayRef to the forces pointer passed in. Virial/energy are also outputs if they were to be computed.

Note: also launches the reinitalization of the PME output buffers. TODO: this should be moved out to avoid miscounting its wall-time (as wait iso launch).

Parameters
[in]pmeThe PME data structure.
[in]wcycleThe wallclock counter.
[out]forcesThe output forces.
[out]virialThe output virial matrix.
[out]energyThe output energy.
[in]completionKindIndicates whether PME task completion should only be checked rather than waited for
Returns
True if the PME GPU tasks have completed
void pme_gpu_wait_finish_task ( const gmx_pme_t *  pme,
gmx_wallcycle_t  wcycle,
gmx::ArrayRef< const gmx::RVec > *  forces,
matrix  virial,
real energy 
)

Blocks until PME GPU tasks are completed, and gets the output forces and virial/energy (if they were to be computed).

Parameters
[in]pmeThe PME data structure.
[out]wcycleThe wallclock counter.
[out]forcesThe output forces.
[out]virialThe output virial matrix.
[out]energyThe output energy.
PmeRunMode pme_run_mode ( const gmx_pme_t *  pme)

Returns the active PME codepath (CPU, GPU, mixed).

Todo:
This is a rather static data that should be managed by the higher level task scheduler.
Parameters
[in]pmeThe PME data structure.
Returns
active PME codepath.