1 #ifndef OPENMM_NONBONDEDFORCE_H_ 
    2 #define OPENMM_NONBONDEDFORCE_H_ 
   41 #include "internal/windowsExport.h" 
   96         CutoffNonPeriodic = 1,
 
  122         return particles.size();
 
  128         return exceptions.size();
 
  133     NonbondedMethod getNonbondedMethod() 
const;
 
  137     void setNonbondedMethod(NonbondedMethod method);
 
  144     double getCutoffDistance() 
const;
 
  151     void setCutoffDistance(
double distance);
 
  156     bool getUseSwitchingFunction() 
const;
 
  161     void setUseSwitchingFunction(
bool use);
 
  166     double getSwitchingDistance() 
const;
 
  171     void setSwitchingDistance(
double distance);
 
  175     double getReactionFieldDielectric() 
const;
 
  179     void setReactionFieldDielectric(
double dielectric);
 
  189     double getEwaldErrorTolerance() 
const;
 
  199     void setEwaldErrorTolerance(
double tol);
 
  209     void getPMEParameters(
double& alpha, 
int& nx, 
int& ny, 
int& nz) 
const;
 
  219     void setPMEParameters(
double alpha, 
int nx, 
int ny, 
int nz);
 
  232     int addParticle(
double charge, 
double sigma, 
double epsilon);
 
  241     void getParticleParameters(
int index, 
double& charge, 
double& sigma, 
double& epsilon) 
const;
 
  252     void setParticleParameters(
int index, 
double charge, 
double sigma, 
double epsilon);
 
  269     int addException(
int particle1, 
int particle2, 
double chargeProd, 
double sigma, 
double epsilon, 
bool replace = 
false);
 
  280     void getExceptionParameters(
int index, 
int& particle1, 
int& particle2, 
double& chargeProd, 
double& sigma, 
double& epsilon) 
const;
 
  293     void setExceptionParameters(
int index, 
int particle1, 
int particle2, 
double chargeProd, 
double sigma, 
double epsilon);
 
  306     void createExceptionsFromBonds(
const std::vector<std::pair<int, int> >& bonds, 
double coulomb14Scale, 
double lj14Scale);
 
  314         return useDispersionCorrection;
 
  323         useDispersionCorrection = useCorrection;
 
  331     int getReciprocalSpaceForceGroup() 
const;
 
  341     void setReciprocalSpaceForceGroup(
int group);
 
  354     void updateParametersInContext(
Context& context);
 
  367     ForceImpl* createImpl() 
const;
 
  371     NonbondedMethod nonbondedMethod;
 
  372     double cutoffDistance, switchingDistance, rfDielectric, ewaldErrorTol, alpha;
 
  373     bool useSwitchingFunction, useDispersionCorrection;
 
  374     int recipForceGroup, nx, ny, nz;
 
  375     void addExclusionsToSet(
const std::vector<std::set<int> >& bonded12, std::set<int>& exclusions, 
int baseParticle, 
int fromParticle, 
int currentLevel) 
const;
 
  376     std::vector<ParticleInfo> particles;
 
  377     std::vector<ExceptionInfo> exceptions;
 
  378     std::map<std::pair<int, int>, 
int> exceptionMap;
 
  385 class NonbondedForce::ParticleInfo {
 
  387     double charge, sigma, epsilon;
 
  389         charge = sigma = epsilon = 0.0;
 
  391     ParticleInfo(
double charge, 
double sigma, 
double epsilon) :
 
  392         charge(charge), sigma(sigma), epsilon(epsilon) {
 
  400 class NonbondedForce::ExceptionInfo {
 
  402     int particle1, particle2;
 
  403     double chargeProd, sigma, epsilon;
 
  405         particle1 = particle2 = -1;
 
  406         chargeProd = sigma = epsilon = 0.0;
 
  408     ExceptionInfo(
int particle1, 
int particle2, 
double chargeProd, 
double sigma, 
double epsilon) :
 
  409         particle1(particle1), particle2(particle2), chargeProd(chargeProd), sigma(sigma), epsilon(epsilon) {
 
A Context stores the complete state of a simulation. 
Definition: Context.h:67
Periodic boundary conditions are used, so that each particle interacts only with the nearest periodic...
Definition: NonbondedForce.h:102
Periodic boundary conditions are used, and Ewald summation is used to compute the interaction of each...
Definition: NonbondedForce.h:107
Periodic boundary conditions are used, and Particle-Mesh Ewald (PME) summation is used to compute the...
Definition: NonbondedForce.h:112
bool getUseDispersionCorrection() const 
Get whether to add a contribution to the energy that approximately represents the effect of Lennard-J...
Definition: NonbondedForce.h:313
Force objects apply forces to the particles in a System, or alter their behavior in other ways...
Definition: Force.h:65
void setUseDispersionCorrection(bool useCorrection)
Set whether to add a contribution to the energy that approximately represents the effect of Lennard-J...
Definition: NonbondedForce.h:322
int getNumExceptions() const 
Get the number of special interactions that should be calculated differently from other interactions...
Definition: NonbondedForce.h:127
int getNumParticles() const 
Get the number of particles for which force field parameters have been defined. 
Definition: NonbondedForce.h:121
bool usesPeriodicBoundaryConditions() const 
Returns whether or not this force makes use of periodic boundary conditions. 
Definition: NonbondedForce.h:361
NonbondedMethod
This is an enumeration of the different methods that may be used for handling long range nonbonded fo...
Definition: NonbondedForce.h:86
This class implements nonbonded interactions between particles, including a Coulomb force to represen...
Definition: NonbondedForce.h:81
Definition: AndersenThermostat.h:40