STORMM Source Documentation
|
A collection of contol data for molecular mechanics simulations, conveying the current step number, progress counters through various work units, the time step, RATTLE tolerance, and other critical parameters to guide calculations. This common struct accommodates both molecular dynamics and molecular mechanics energy minimizations, and can be created from control objects derived from &dynamics or &minimize namelists. More...
#include <mm_controls.h>
Public Member Functions | |
MolecularMechanicsControls (const MolecularMechanicsControls &original) | |
The copy constructor handles assignment of internal POINTER-kind Hybrid objects. | |
MolecularMechanicsControls (MolecularMechanicsControls &&original) | |
The move constructor prepared the original object for destruction. | |
MolecularMechanicsControls & | operator= (const MolecularMechanicsControls &other) |
The copy assignment operator likewise handles assignment of internal POINTER-kind Hybrid objects. | |
MolecularMechanicsControls & | operator= (MolecularMechanicsControls &&other) |
The move assignment operator looks much like the copy assignment operator. | |
int | getStepNumber () const |
Get the current step number. | |
int | getSteepestDescentCycles () const |
Get the number of steepest descent minimization cycles. | |
int | getTotalCycles () const |
Get the total number of minimization steps, or total MD cycles. | |
int | getNTWarpMultiplicity () const |
Get the neutral-territory warp multiplicity. | |
double | getInitialMinimizationStep () const |
Get the initial step for energy minimization. | |
double | getElectrostaticCutoff () const |
Get the cutoff for non-bonded electrostatic interactions in periodic simulations. | |
double | getVanDerWaalsCutoff () const |
Get the cutoff for non-bonded van-der Waals interactions in periodic simulations. | |
int | getValenceWorkUnitProgress (int counter_index, HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get the value of one of the valence work unit progress counters on the host or the HPC device. | |
int | getNonbondedWorkUnitProgress (int counter_index, HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get the value of one of the non-bonded work unit progress counters on the host or the HPC device. | |
int | getPmeWorkUnitProgress (int counter_index, HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get the value of one of the Particle-Mesh Ewald work unit progress counters on the host or the HPC device. | |
int | getReductionWorkUnitProgress (int counter_index, ReductionStage process, HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get the value of one of the reduction work unit progress counters on the host or the HPC device. | |
MMControlKit< double > | dpData (HybridTargetLevel tier=HybridTargetLevel::HOST) |
Obtain a double-precision abstract for this object. | |
MMControlKit< float > | spData (HybridTargetLevel tier=HybridTargetLevel::HOST) |
Obtain a single-precision abstract for this object. | |
void | incrementStep () |
Increment the step counter, moving the controls to a different progress counter. | |
MolecularMechanicsControls (double initial_step_in=default_minimize_dx0, int sd_cycles_in=default_minimize_ncyc, int max_cycles_in=default_minimize_maxcyc, int nt_warp_multiplicity_in=default_nt_warp_multiplicity, double electrostatic_cutoff_in=default_electrostatic_cutoff, double van_der_waals_cutoff_in=default_van_der_waals_cutoff) | |
The constructor can create an empty object with default parameters for the time step and rattle tolerance, accept user-specified values for critical constants, or be built from a namelist-derived object on the CPU. The reason that those namelist-derived control objects don't just get pushed by value to the GPU is that the work unit counters need special arrays, which this object manages. | |
MolecularMechanicsControls (const DynamicsControls &user_input) | |
MolecularMechanicsControls (const MinimizeControls &user_input) | |
void | primeWorkUnitCounters (const CoreKlManager &launcher, EvaluateForce eval_frc, EvaluateEnergy eval_nrg, ClashResponse softcore, VwuGoal purpose, PrecisionModel valence_prec, PrecisionModel nonbond_prec, QMapMethod qspread_approach, PrecisionModel acc_prec, size_t image_coord_type, int qspread_order, NeighborListKind nbgr_config, TinyBoxPresence has_tiny_box, const AtomGraphSynthesis &poly_ag) |
Prime the work unit counters based on a particular GPU configuration. | |
void | primeWorkUnitCounters (const CoreKlManager &launcher, EvaluateForce eval_frc, EvaluateEnergy eval_nrg, const ClashResponse softcore, VwuGoal purpose, PrecisionModel valence_prec, PrecisionModel nonbond_prec, const AtomGraphSynthesis &poly_ag) |
void | primeWorkUnitCounters (const CoreKlManager &launcher, EvaluateForce eval_frc, EvaluateEnergy eval_nrg, VwuGoal purpose, PrecisionModel valence_prec, PrecisionModel nonbond_prec, const AtomGraphSynthesis &poly_ag) |
void | primeWorkUnitCounters (const CoreKlManager &launcher, EvaluateForce eval_frc, EvaluateEnergy eval_nrg, VwuGoal purpose, PrecisionModel general_prec, const AtomGraphSynthesis &poly_ag) |
template<typename T, typename Tacc, typename Tcalc, typename T4> | |
void | setNTWarpMultiplicity (const CellGrid< T, Tacc, Tcalc, T4 > *cg_a, const CellGrid< T, Tacc, Tcalc, T4 > *cg_b, const GpuDetails &gpu) |
Set the neutral-territory warp multiplicity based on one or two gell grids, for a given GPU. | |
template<typename T, typename Tacc, typename Tcalc, typename T4> | |
void | setNTWarpMultiplicity (const CellGrid< T, Tacc, Tcalc, T4 > *cg_a, const GpuDetails &gpu) |
void | setNTWarpMultiplicity (int mult_in) |
A collection of contol data for molecular mechanics simulations, conveying the current step number, progress counters through various work units, the time step, RATTLE tolerance, and other critical parameters to guide calculations. This common struct accommodates both molecular dynamics and molecular mechanics energy minimizations, and can be created from control objects derived from &dynamics or &minimize namelists.
stormm::mm::MolecularMechanicsControls::MolecularMechanicsControls | ( | double | initial_step_in = default_minimize_dx0, |
int | sd_cycles_in = default_minimize_ncyc, | ||
int | max_cycles_in = default_minimize_maxcyc, | ||
int | nt_warp_multiplicity_in = default_nt_warp_multiplicity, | ||
double | electrostatic_cutoff_in = default_electrostatic_cutoff, | ||
double | van_der_waals_cutoff_in = default_van_der_waals_cutoff ) |
The constructor can create an empty object with default parameters for the time step and rattle tolerance, accept user-specified values for critical constants, or be built from a namelist-derived object on the CPU. The reason that those namelist-derived control objects don't just get pushed by value to the GPU is that the work unit counters need special arrays, which this object manages.
initial_step_in | The desired initial step (for energy minimization) |
user_input | Namelist-derived molecular dynamics or energy minimization controls |
stormm::mm::MolecularMechanicsControls::MolecularMechanicsControls | ( | const MolecularMechanicsControls & | original | ) |
The copy constructor handles assignment of internal POINTER-kind Hybrid objects.
original | The object to copy |
stormm::mm::MolecularMechanicsControls::MolecularMechanicsControls | ( | MolecularMechanicsControls && | original | ) |
The move constructor prepared the original object for destruction.
original | The object from which to take contents |
int stormm::mm::MolecularMechanicsControls::getNonbondedWorkUnitProgress | ( | int | counter_index, |
HybridTargetLevel | tier = HybridTargetLevel::HOST ) const |
Get the value of one of the non-bonded work unit progress counters on the host or the HPC device.
counter_index | Index of the work unit counter to retrieve |
tier | Get the result from the CPU host or the HPC device |
int stormm::mm::MolecularMechanicsControls::getPmeWorkUnitProgress | ( | int | counter_index, |
HybridTargetLevel | tier = HybridTargetLevel::HOST ) const |
Get the value of one of the Particle-Mesh Ewald work unit progress counters on the host or the HPC device.
counter_index | Index of the work unit counter to retrieve |
tier | Get the result from the CPU host or the HPC device |
int stormm::mm::MolecularMechanicsControls::getReductionWorkUnitProgress | ( | int | counter_index, |
ReductionStage | process, | ||
HybridTargetLevel | tier = HybridTargetLevel::HOST ) const |
Get the value of one of the reduction work unit progress counters on the host or the HPC device.
counter_index | Index of the work unit counter to retrieve |
tier | Get the result from the CPU host or the HPC device |
int stormm::mm::MolecularMechanicsControls::getValenceWorkUnitProgress | ( | int | counter_index, |
HybridTargetLevel | tier = HybridTargetLevel::HOST ) const |
Get the value of one of the valence work unit progress counters on the host or the HPC device.
counter_index | Index of the work unit counter to retrieve |
tier | Get the result from the CPU host or the HPC device |
MolecularMechanicsControls & stormm::mm::MolecularMechanicsControls::operator= | ( | const MolecularMechanicsControls & | other | ) |
The copy assignment operator likewise handles assignment of internal POINTER-kind Hybrid objects.
other | Another way to say original, in a different semantic context |
MolecularMechanicsControls & stormm::mm::MolecularMechanicsControls::operator= | ( | MolecularMechanicsControls && | other | ) |
The move assignment operator looks much like the copy assignment operator.
other | The object from which to take contents |
void stormm::mm::MolecularMechanicsControls::primeWorkUnitCounters | ( | const CoreKlManager & | launcher, |
EvaluateForce | eval_frc, | ||
EvaluateEnergy | eval_nrg, | ||
ClashResponse | softcore, | ||
VwuGoal | purpose, | ||
PrecisionModel | valence_prec, | ||
PrecisionModel | nonbond_prec, | ||
QMapMethod | qspread_approach, | ||
PrecisionModel | acc_prec, | ||
size_t | image_coord_type, | ||
int | qspread_order, | ||
NeighborListKind | nbgr_config, | ||
TinyBoxPresence | has_tiny_box, | ||
const AtomGraphSynthesis & | poly_ag ) |
Prime the work unit counters based on a particular GPU configuration.
Overloaded:
launcher | Object containing launch parameters for all kernels |
eval_frc | Indicate whether forces are to be evaluated–some kernels allocate different numbers of blocks in the launch grid if forces are not required. |
eval_nrg | Indicate whether energy is to be evaluated–some kernels allocate different numbers of blocks in the launch grid if the energy is not required. |
softcore | Indicate whether the control counters should take clash handling kernel variants into account |
valence_prec | Precision model for valence calculations |
nonbond_prec | Precision model for non-bonded calculations |
general_prec | Precision model for all calculations |
qspread_approach | The approach taken by GPU kernels to spread charge density onto the particle-mesh interaction grids |
acc_prec | Precision model for accumulations (needed in case the density mapping approach is accumulation in the shared memory space) |
image_coord_type | Data type used to represent coordinates in the neighbor list images |
qspread_order | The order of B-spline interpolation used to spread charge density to the particle-mesh interaction grids |
poly_ag | Compilation of topologies describing the workload (used here for general descriptors such as the non-bonded work unit type) |
void stormm::mm::MolecularMechanicsControls::setNTWarpMultiplicity | ( | const CellGrid< T, Tacc, Tcalc, T4 > * | cg_a, |
const CellGrid< T, Tacc, Tcalc, T4 > * | cg_b, | ||
const GpuDetails & | gpu ) |
Set the neutral-territory warp multiplicity based on one or two gell grids, for a given GPU.
Overloaded:
cg_a | The first cell grid neighbor list |
cg_b | The second (optional) cell grid neighbor list |
gpu | Details of the GPU that will evaluate pairwise interactions |
mult_in | The multiplicity to explicitly set |