Highly Optimized Object-oriented Many-particle Dynamics -- Blue Edition

Lennard-Jones pair force. More...

Public Member Functions

def __init__
 Specify the Lennard-Jones pair force. More...
 
def set_params
 Set parameters controlling the way forces are computed. More...
 
def disable
 Disables the force. More...
 
def benchmark
 Benchmarks the force computation. More...
 
def enable
 Enables the force. More...
 

Detailed Description

Lennard-Jones pair force.

The command pair.lj specifies that a Lennard-Jones type pair force should be added to every non-bonded particle pair in the simulation.

\begin{eqnarray*} V_{\mathrm{LJ}}(r) = & 4 \varepsilon \left[ \left( \frac{\sigma}{r} \right)^{12} - \alpha \left( \frac{\sigma}{r} \right)^{6} \right] & r < r_{\mathrm{cut}} \\ = & 0 & r \ge r_{\mathrm{cut}} \\ \end{eqnarray*}

For an exact definition of the force and potential calculation and how cutoff radii are handled, see pair.

The following coefficients must be set per unique pair of particle types. See hoomd_script.pair or the Quick Start Tutorial for information on how to set coefficients.

  • \( \varepsilon \) - epsilon (in energy units)
  • \( \sigma \) - sigma (in distance units)
  • \( \alpha \) - alpha (unitless)
    • optional: defaults to 1.0
  • \( r_{\mathrm{cut}} \) - r_cut (in distance units)
    • optional: defaults to the global r_cut specified in the pair command
  • \( r_{\mathrm{on}} \) - r_on (in distance units)
    • optional: defaults to the global r_cut specified in the pair command

pair.lj is a standard pair potential and supports a number of energy shift / smoothing modes. See hoomd_script.pair.pair for a full description of the various options.

Example:

1 lj.pair_coeff.set('A', 'A', epsilon=1.0, sigma=1.0)
2 lj.pair_coeff.set('A', 'B', epsilon=2.0, sigma=1.0, alpha=0.5, r_cut=3.0, r_on=2.0);
3 lj.pair_coeff.set('B', 'B', epsilon=1.0, sigma=1.0, r_cut=2**(1.0/6.0), r_on=2.0);
4 lj.pair_coeff.set(['A', 'B'], ['C', 'D'], epsilon=1.5, sigma=2.0)

For more information on setting pair coefficients, including examples with wildcards, see pair_coeff.set().

The cutoff radius r_cut passed into the initial pair.lj command sets the default r_cut for all pair interactions. Smaller (or larger) cutoffs can be set individually per each type pair. The cutoff distances used for the neighbor list will by dynamically determined from the maximum of all r_cut values specified among all type pair parameters among all pair potentials.

MPI Support
This command works in MPI parallel jobs.
Examples:
analyze_imd, create_random_polymers, dump_dcd, init_create_empty, init_reset, and init_xml.

Constructor & Destructor Documentation

def __init__ (   self,
  r_cut,
  name = None 
)

Specify the Lennard-Jones pair force.

Parameters
r_cutDefault cutoff radius (in distance units)
nameName of the force instance

Example:

1 lj = pair.lj(r_cut=3.0)
2 lj.pair_coeff.set('A', 'A', epsilon=1.0, sigma=1.0)
3 lj.pair_coeff.set('A', 'B', epsilon=2.0, sigma=1.0, alpha=0.5, r_cut=3.0, r_on=2.0);
4 lj.pair_coeff.set('B', 'B', epsilon=1.0, sigma=1.0, r_cut=2**(1.0/6.0), r_on=2.0);
Note
Pair coefficients for all type pairs in the simulation must be set before it can be started with run()

Member Function Documentation

def benchmark (   self,
  n 
)
inherited

Benchmarks the force computation.

Parameters
nNumber of iterations to average the benchmark over

Examples:

1 t = force.benchmark(n = 100)

The value returned by benchmark() is the average time to perform the force computation, in milliseconds. The benchmark is performed by taking the current positions of all particles in the simulation and repeatedly calculating the forces on them. Thus, you can benchmark different situations as you need to by simply running a simulation to achieve the desired state before running benchmark().

Note
There is, however, one subtle side effect. If the benchmark() command is run directly after the particle data is initialized with an init command, then the results of the benchmark will not be typical of the time needed during the actual simulation. Particles are not reordered to improve cache performance until at least one time step is performed. Executing run(1) before the benchmark will solve this problem.

To use this command, you must have saved the force in a variable, as shown in this example:

1 force = pair.some_force()
2 # ... later in the script
3 t = force.benchmark(n = 100)
def disable (   self,
  log = False 
)
inherited

Disables the force.

Parameters
logSet to True if you plan to continue logging the potential energy associated with this force.

Examples:

1 force.disable()
2 force.disable(log=True)

Executing the disable command will remove the force from the simulation. Any run() command executed after disabling a force will not calculate or use the force during the simulation. A disabled force can be re-enabled with enable()

By setting log to True, the values of the force can be logged even though the forces are not applied in the simulation. For forces that use cutoff radii, setting log=True will cause the correct r_cut values to be used throughout the simulation, and therefore possibly drive the neighbor list size larger than it otherwise would be. If log is left False, the potential energy associated with this force will not be available for logging.

To use this command, you must have saved the force in a variable, as shown in this example:

1 force = pair.some_force()
2 # ... later in the script
3 force.disable()
4 force.disable(log=True)
def enable (   self)
inherited

Enables the force.

Examples:

1 force.enable()

See disable() for a detailed description.

def set_params (   self,
  mode = None 
)
inherited

Set parameters controlling the way forces are computed.

Parameters
mode(if set) Set the mode with which potentials are handled at the cutoff

valid values for mode are: "none" (the default), "shift", and "xplor"

  • none - No shifting is performed and potentials are abruptly cut off
  • shift - A constant shift is applied to the entire potential so that it is 0 at the cutoff
  • xplor - A smoothing function is applied to gradually decrease both the force and potential to 0 at the cutoff when ron < rcut, and shifts the potential to 0 ar the cutoff when ron >= rcut. (see pair above for formulas and more information)

Examples:

1 mypair.set_params(mode="shift")
2 mypair.set_params(mode="no_shift")
3 mypair.set_params(mode="xplor")

Documentation for HOOMD-blue v0.11.3-1378-ge1e0064
Generated on Fri Apr 18 2014 07:30:07 for HOOMD-blue by doxygen 1.8.5