STORMM Source Documentation
Loading...
Searching...
No Matches
stormm::energy::PMIGrid Class Reference

An object to hold a series of meshes for accumulating density from condensed-phase molecular systems and transforming it into potential, whether for charges and electrostatics or dispersion sources and van-der Waals interactions. This supports the particle-mesh interactions (PMI) of the particle-particle / particle-mesh potential partitioning common in several simulation techniques. More...

#include <pmigrid.h>

Public Member Functions

NonbondedTheme getTheme () const
 Get the non-bonded propery mapped onto the particle-mesh interaction grid.
 
PrecisionModel getMode () const
 Get the mode in which the object is set to operate.
 
FFTMode getFFTStaging () const
 Get the method for performing real-to-complex and complex-to-real FFTs. This implies how the grid data itself may be padded.
 
bool fixedPrecisionEnabled () const
 Indicate whether fixed-precision accumulation is enabled. This will test whether fp_accumulation_bits is greater than zero.
 
QMapMethod getRecommendedMappingMethod () const
 Get the recommended density mapping method for the current set of grids.
 
QMapMethod getWorkUnitConfiguration () const
 Get the work unit configuration of the object.
 
bool dataIsReal () const
 Indicate whether the data is present as real numbers or fixed-precision accumulators. Even with fixed-precision enabled, this object may hold data as real numbers for non-accumulating tasks.
 
int getFixedPrecisionBits () const
 Get the number of bits used in scaling numbers for fixed precision accumulation.
 
int getSharedFixedPrecisionBits () const
 Get the number of bits used in scaling numbers for fixed precision accumulation in shared accumulation density mapping kernels.
 
int getSystemCount () const
 Get the number of systems with grid managed by this object.
 
int getInterpolationOrder () const
 Get the interpolation order to be used in mapping particle density to the particle-mesh interaction grid.
 
size_t getTotalCapacity () const
 Get the total capacity allocated to store grid data.
 
int getWorkUnitCount () const
 Get the number of work units used by the optimized kernel.
 
int getLargestWorkUnitGridPoints () const
 Report the number of grid points present in the largest work unit.
 
bool shortFormatAccumulation () const
 Report whether the PMIGrid will use short format accumulation.
 
bool useOverflowAccumulation () const
 Report whether the PMIGrid will need to use overflow accumulators.
 
template<typename T, typename Tacc, typename Tcalc, typename T4>
const CellGrid< T, Tacc, Tcalc, T4 > * getCellGridPointer () const
 Get the pointer to the attached CellGrid object, which provides the basis for the grids / meshes contained in this object. This will restore the type of the object but must be invoked in a manner consistent with its original typing.
 
const CellGridReader< void, void, void, void > getTemplateFreeCellGridReader (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get an abstract of the attached cell grid object, with all templating cast away. This exists because the templated types of the pointer stored in the PMIGrid object are placeholders for its true types. Out of an abundance of caution, a CellGrid pointer with the proper data types for all four components is created so that the CellGrid's public member function templateFreeData() can be invoked.
 
size_t getCellGridMatrixTypeID () const
 Get the type ID of the attached CellGrid object's cell dimension matrix expressions. This is the "typename T" seen in many CellGrid class function declarations, and the root data type of the four-tuples in which its coordinates are expressed. Possible types include float, double, int, and llint.
 
size_t getCellGridAccumulatorTypeID () const
 Get the type ID of the attached CellGrid object's accumulators. This is "typename Tacc" in various CellGrid class function declarations and invocations. Acceptable types include int and llint.
 
size_t getCellGridCalculationTypeID () const
 Get the type ID of the attached CellGrid object's fractional space coordinate transformations, which sets a precedent for any calculations that will be done with an object of the class. This is "typename Tcalc" in various CellGrid class function declarations and invocations.
 
size_t getCellGridCoordinateTypeID () const
 Get the type ID of the attached CellGrid object's coordinate data, a four-tuple of its cell dimension matrix data. Acceptable types include float, double, int, and llint. This is "typename T4" in most CellGrid-associated function declarations and invocations, sometimes "typename Tcrd".
 
double getTotalOnGrid (int system_index) const
 Get the total substance, particle density, mapped to all grid elements for a particular system.
 
std::vector< double > getGrid (int system_index) const
 Export the particle-mesh interaction grid asssociated with one system in the synthesis. The resulting linearized data reflects an ordering of the grid that is identical to the way it appears in the PMIGrid itself.
 
const PhaseSpaceSynthesisgetCoordinateSynthesisPointer () const
 Get a pointer to the underlying coordinate synthesis.
 
const AtomGraphSynthesisgetTopologySynthesisPointer () const
 Get a pointer to the underlying topology synthesis.
 
const PMIGridgetSelfPointer () const
 Get a pointer to the PMIGrid object itself.
 
void setMode (PrecisionModel mode_in)
 Change the mode in which the object is set to operate. This will leave bounds arrays unchanged but blow away data as the object is reconfigured.
 
void setRealDataFormat (bool real_in=true)
 Change the character of the data format. If called with TRUE the object will report holding real-valued data.
 
void prepareFixedPrecisionModel (int fp_accumulation_bits_in, int shared_fp_bits_in)
 Set the object to enable or disable fixed-precision accumulation.
 
void setRecommendedMappingMethod (const GpuDetails &gpu)
 Set the recommended mapping method based on the density of the grid and available resources.
 
void initialize (HybridTargetLevel tier=HybridTargetLevel::HOST, const GpuDetails &gpu=null_gpu)
 Initialize the particle-mesh interaction grids, respecting the object's precision and accumulation modes. If fixed-precision is enabled, this will set the data content to register as no longer "real numbers" (not imaginary numbers, but integers).
 
void convertToReal (HybridTargetLevel tier=HybridTargetLevel::HOST, const GpuDetails &gpu=null_gpu)
 Convert fixed-precision integer data in the object to real data. There is no function to convert real back to fixed-precision integer data: the process is to accumulate in the one mode and then convert to the other when finished.
 
template<typename T, typename Tacc, typename Tcalc, typename T4>
 PMIGrid (const CellGrid< T, Tacc, Tcalc, T4 > *cg_in, NonbondedTheme theme_in, int b_spline_order_in=default_bspline_order, PrecisionModel mode_in=PrecisionModel::SINGLE, FFTMode fft_staging_in=FFTMode::IN_PLACE, int fp_accumulation_bits_in=0, int shared_fp_accumulation_bits_in=-1, const GpuDetails &gpu=null_gpu, const QMapMethod work_unit_configuration_in=QMapMethod::GENERAL_PURPOSE)
 The constructor builds the object off of a CellGrid object, which is in turn linked to a particular coordinate synthesis (PhaseSpaceSynthesis).
 
template<typename T, typename Tacc, typename Tcalc, typename T4>
 PMIGrid (const CellGrid< T, Tacc, Tcalc, T4 > &cg_in, NonbondedTheme theme_in, int b_spline_order_in=default_bspline_order, PrecisionModel mode_in=PrecisionModel::SINGLE, FFTMode fft_staging_in=FFTMode::IN_PLACE, int fp_accumulation_bits_in=0, int shared_fp_accumulation_bits_in=-1, const GpuDetails &gpu=null_gpu, const QMapMethod work_unit_configuration_in=QMapMethod::GENERAL_PURPOSE)
 
 PMIGrid (const PMIGrid &original)=default
 The choice to re-cast arrays in the abstracts rather than have special POINTER-kind Hybrid objects in the PMIGrid itself allows the default copy and move constructors, as well as the copy and move assignment operators, to remain valid.
 
 PMIGrid (PMIGrid &&original)=default
 
PMIGridoperator= (const PMIGrid &original)=default
 
PMIGridoperator= (PMIGrid &&original)=default
 
uint4 getGridDimensions (int system_index) const
 Obtain the dimensions of one of the grids in the object.
 
int getGridDimensions (int system_index, UnitCellAxis uc_axis) const
 
int getGridDimensions (int system_index, CartesianDimension uc_axis) const
 
PMIGridWriter data (HybridTargetLevel tier=HybridTargetLevel::HOST)
 Get the abstract for this object containing its precision mode, grid dimensions, and pointers to relevant data arrays.
 
const PMIGridReader data (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
PMIGridAccumulator fpData (HybridTargetLevel tier=HybridTargetLevel::HOST)
 Get the fixed-precision abstract for this object, re-interpreting either main data array and adding an integer array of overflow accumulators for accumulation protected against the non-associativity that breaks bitwise reproducibility for floating point parallel computations.
 
const PMIGridFPReader fpData (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
void prepareWorkUnits (QMapMethod approach, const GpuDetails &gpu)
 Prepare a set of work units for the synthesis of systems. These work units will be tailored to the synthesis for which the topology was created, but can be changed given a new GPU.
 
void prepareWorkUnits (const GpuDetails &gpu)
 

Detailed Description

An object to hold a series of meshes for accumulating density from condensed-phase molecular systems and transforming it into potential, whether for charges and electrostatics or dispersion sources and van-der Waals interactions. This supports the particle-mesh interactions (PMI) of the particle-particle / particle-mesh potential partitioning common in several simulation techniques.

Constructor & Destructor Documentation

◆ PMIGrid()

stormm::energy::PMIGrid::PMIGrid ( const PMIGrid & original)
default

The choice to re-cast arrays in the abstracts rather than have special POINTER-kind Hybrid objects in the PMIGrid itself allows the default copy and move constructors, as well as the copy and move assignment operators, to remain valid.

Parameters
originalThe original object to copy or move
otherAnother object placed on the right hand side of the assignment statement

Member Function Documentation

◆ convertToReal()

void stormm::energy::PMIGrid::convertToReal ( HybridTargetLevel tier = HybridTargetLevel::HOST,
const GpuDetails & gpu = null_gpu )

Convert fixed-precision integer data in the object to real data. There is no function to convert real back to fixed-precision integer data: the process is to accumulate in the one mode and then convert to the other when finished.

Parameters
tierIndicate whether to initialize the density grids at the level of the CPU host or the GPU device
gpuDetails of the available GPU, if the initialization focuses on device memory

◆ data()

PMIGridWriter stormm::energy::PMIGrid::data ( HybridTargetLevel tier = HybridTargetLevel::HOST)

Get the abstract for this object containing its precision mode, grid dimensions, and pointers to relevant data arrays.

Overloaded:

  • Get a writeable abstract from a mutable object.
  • Get a read-only abstract from a const object.
Parameters
tierSpecify whether to obtain pointers on the CPU host or GPU device

◆ fpData()

PMIGridAccumulator stormm::energy::PMIGrid::fpData ( HybridTargetLevel tier = HybridTargetLevel::HOST)

Get the fixed-precision abstract for this object, re-interpreting either main data array and adding an integer array of overflow accumulators for accumulation protected against the non-associativity that breaks bitwise reproducibility for floating point parallel computations.

Overloaded:

  • Get a writeable abstract from a mutable object.
  • Get a read-only abstract from a const object.
Parameters
tierSpecify whether to obtain pointers on the CPU host or GPU device

◆ getGrid()

std::vector< double > stormm::energy::PMIGrid::getGrid ( int system_index) const

Export the particle-mesh interaction grid asssociated with one system in the synthesis. The resulting linearized data reflects an ordering of the grid that is identical to the way it appears in the PMIGrid itself.

Parameters
system_indexThe system of interest

◆ getGridDimensions()

uint4 stormm::energy::PMIGrid::getGridDimensions ( int system_index) const

Obtain the dimensions of one of the grids in the object.

Overloaded:

  • Get all dimensions, including the offset in the data arrays, of one system's PMI grid
  • Get the dimension of one system's grid along a particular axis
Parameters
system_indexIndex of the system of interest within the associated synthesis
uc_axisThe unit cell axis of interest

◆ getTemplateFreeCellGridReader()

const CellGridReader< void, void, void, void > stormm::energy::PMIGrid::getTemplateFreeCellGridReader ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get an abstract of the attached cell grid object, with all templating cast away. This exists because the templated types of the pointer stored in the PMIGrid object are placeholders for its true types. Out of an abundance of caution, a CellGrid pointer with the proper data types for all four components is created so that the CellGrid's public member function templateFreeData() can be invoked.

Parameters
tierSpecify whether to obtain pointers, void-casted or otherwise, on the CPU host or GPU device

◆ getTotalOnGrid()

double stormm::energy::PMIGrid::getTotalOnGrid ( int system_index) const

Get the total substance, particle density, mapped to all grid elements for a particular system.

Parameters
system_indexThe system of interest

◆ initialize()

void stormm::energy::PMIGrid::initialize ( HybridTargetLevel tier = HybridTargetLevel::HOST,
const GpuDetails & gpu = null_gpu )

Initialize the particle-mesh interaction grids, respecting the object's precision and accumulation modes. If fixed-precision is enabled, this will set the data content to register as no longer "real numbers" (not imaginary numbers, but integers).

Parameters
tierIndicate whether to initialize the density grids at the level of the CPU host or the GPU device
gpuDetails of the available GPU, if the initialization focuses on device memory

◆ prepareFixedPrecisionModel()

void stormm::energy::PMIGrid::prepareFixedPrecisionModel ( int fp_accumulation_bits_in,
int shared_fp_bits_in )

Set the object to enable or disable fixed-precision accumulation.

Parameters
fp_accumulation_bits_inThe number of bits used in fixed-precision accumulation, indicating the power of two by which to scale numbers going int each sum
shared_fp_bits_inThe number of bits used in fixed-precision accumulation with a separate space in shared memory. While this will be set to be consistent with fp_accumulation_bits_in if possible, there is a critical difference: this accumulation is done in the interest of going directly to

◆ prepareWorkUnits()

void stormm::energy::PMIGrid::prepareWorkUnits ( QMapMethod approach,
const GpuDetails & gpu )

Prepare a set of work units for the synthesis of systems. These work units will be tailored to the synthesis for which the topology was created, but can be changed given a new GPU.

Overloaded:

  • Tailor work units to a chosen a mapping method
  • Tailor work units to the recommended mapping method
Parameters
approachThe chosen work unit mapping method (if the grid density or parameters do do not comport with the method, an alternative will be selected)
gpuDetails of the available GPU

◆ setMode()

void stormm::energy::PMIGrid::setMode ( PrecisionModel mode_in)

Change the mode in which the object is set to operate. This will leave bounds arrays unchanged but blow away data as the object is reconfigured.

Parameters
mode_inThe chosen mode for the object

◆ setRealDataFormat()

void stormm::energy::PMIGrid::setRealDataFormat ( bool real_in = true)

Change the character of the data format. If called with TRUE the object will report holding real-valued data.

Parameters
real_inSpecify whether the object holds real (TRUE) or fixed-precision data (FALSE)

◆ setRecommendedMappingMethod()

void stormm::energy::PMIGrid::setRecommendedMappingMethod ( const GpuDetails & gpu)

Set the recommended mapping method based on the density of the grid and available resources.

Parameters
gpuSpecifications of the GPU that will perform the work

The documentation for this class was generated from the following files: