Update MTK documentation
[hoomd-blue.git] / python-module / hoomd_script / wall.py
blobcc7006aa61ec50dad6489e954a5d1cc23d1cb868
1 # -- start license --
2 # Highly Optimized Object-oriented Many-particle Dynamics -- Blue Edition
3 # (HOOMD-blue) Open Source Software License Copyright 2009-2014 The Regents of
4 # the University of Michigan All rights reserved.
6 # HOOMD-blue may contain modifications ("Contributions") provided, and to which
7 # copyright is held, by various Contributors who have granted The Regents of the
8 # University of Michigan the right to modify and/or distribute such Contributions.
10 # You may redistribute, use, and create derivate works of HOOMD-blue, in source
11 # and binary forms, provided you abide by the following conditions:
13 # * Redistributions of source code must retain the above copyright notice, this
14 # list of conditions, and the following disclaimer both in the code and
15 # prominently in any materials provided with the distribution.
17 # * Redistributions in binary form must reproduce the above copyright notice, this
18 # list of conditions, and the following disclaimer in the documentation and/or
19 # other materials provided with the distribution.
21 # * All publications and presentations based on HOOMD-blue, including any reports
22 # or published results obtained, in whole or in part, with HOOMD-blue, will
23 # acknowledge its use according to the terms posted at the time of submission on:
24 # http://codeblue.umich.edu/hoomd-blue/citations.html
26 # * Any electronic documents citing HOOMD-Blue will link to the HOOMD-Blue website:
27 # http://codeblue.umich.edu/hoomd-blue/
29 # * Apart from the above required attributions, neither the name of the copyright
30 # holder nor the names of HOOMD-blue's contributors may be used to endorse or
31 # promote products derived from this software without specific prior written
32 # permission.
34 # Disclaimer
36 # THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDER AND CONTRIBUTORS ``AS IS'' AND
37 # ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
38 # WARRANTIES OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE, AND/OR ANY
39 # WARRANTIES THAT THIS SOFTWARE IS FREE OF INFRINGEMENT ARE DISCLAIMED.
41 # IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT,
42 # INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
43 # BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
44 # DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
45 # LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE
46 # OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
47 # ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
48 # -- end license --
50 # Maintainer: joaander / All Developers are free to add commands for new features
52 from hoomd_script import force;
53 from hoomd_script import globals;
54 import hoomd;
55 import sys;
56 import math;
57 from hoomd_script import util;
59 ## \package hoomd_script.wall
60 # \brief Commands that specify %wall forces
62 # Walls can add forces to any particles within a certain distance of the wall. Walls are created
63 # when an input XML file is read (read.xml).
65 # By themselves, walls that have been specified in an input file do nothing. Only when you
66 # specify a wall force (i.e. wall.lj), are forces actually applied between the wall and the
67 # particle.
69 ## Lennard-Jones %wall %force
71 # The command wall.lj specifies that a Lennard-Jones type %wall %force should be added to every
72 # particle in the simulation.
74 # The %force \f$ \vec{F}\f$ is
75 # \f{eqnarray*}
76 # \vec{F} = & -\nabla V(r) & r < r_{\mathrm{cut}} \\
77 # = & 0 & r \ge r_{\mathrm{cut}} \\
78 # \f}
79 # where
80 # \f[ V(r) = 4 \varepsilon \left[ \left( \frac{\sigma}{r} \right)^{12} -
81 # \alpha \left( \frac{\sigma}{r} \right)^{6} \right] \f]
82 # and \f$ \vec{r} \f$ is the vector pointing from the %wall to the particle parallel to the wall's normal.
84 # The following coefficients must be set for each particle type using set_coeff().
85 # - \f$ \varepsilon \f$ - \c epsilon (in energy units)
86 # - \f$ \sigma \f$ - \c sigma (in distance units)
87 # - \f$ \alpha \f$ - \c alpha (unitless)
89 # \b Example:
90 # \code
91 # lj.set_coeff('A', epsilon=1.0, sigma=1.0, alpha=1.0)
92 # \endcode
94 # This interaction is applied between every particle and every wall defined in the simulation box. Walls are specified
95 # specified in file given to init.read_xml(). See the page \ref page_xml_file_format for information on creating walls,
96 # specifically the section \ref sec_xml_wall.
98 # The cutoff radius \f$ r_{\mathrm{cut}} \f$ is set once when wall.lj is specified (see __init__())
100 # \MPI_NOT_SUPPORTED
101 class lj(force._force):
102 ## Specify the Lennard-Jones %wall %force
104 # \param r_cut Cutoff radius (in distance units)
106 # \b Example:
107 # \code
108 # lj_wall = wall.lj(r_cut=3.0);
109 # \endcode
111 # \note Coefficients must be set with set_coeff() before the simulation can be run().
112 def __init__(self, r_cut):
113 util.print_status_line();
115 # Error out in MPI simulations
116 if (hoomd.is_MPI_available()):
117 if globals.system_definition.getParticleData().getDomainDecomposition():
118 globals.msg.error("wall.lj is not supported in multi-processor simulations.\n\n")
119 raise RuntimeError("Error setting up wall potential.")
121 # initialize the base class
122 force._force.__init__(self);
124 # create the c++ mirror class
125 self.cpp_force = hoomd.LJWallForceCompute(globals.system_definition, r_cut);
127 # variable for tracking which particle type coefficients have been set
128 self.particle_types_set = [];
130 globals.system.addCompute(self.cpp_force, self.force_name);
132 ## Sets the particle-wall interaction coefficients for a particular particle type
134 # \param particle_type Particle type to set coefficients for
135 # \param epsilon Coefficient \f$ \varepsilon \f$ in the %force (in energy units)
136 # \param sigma Coefficient \f$ \sigma \f$ in the %force (in distance units)
137 # \param alpha Coefficient \f$ \alpha \f$ in the %force (unitless)
139 # Using set_coeff() requires that the specified %wall %force has been saved in a variable. i.e.
140 # \code
141 # lj_wall = wall.lj(r_cut=3.0)
142 # \endcode
144 # \b Examples:
145 # \code
146 # lj_wall.set_coeff('A', epsilon=1.0, sigma=1.0, alpha=1.0)
147 # lj_wall.set_coeff('B', epsilon=1.0, sigma=2.0, alpha=0.0)
148 # \endcode
150 # The coefficients for every particle type in the simulation must be set
151 # before the run() can be started.
152 def set_coeff(self, particle_type, epsilon, sigma, alpha):
153 util.print_status_line();
155 # calculate the parameters
156 lj1 = 4.0 * epsilon * math.pow(sigma, 12.0);
157 lj2 = alpha * 4.0 * epsilon * math.pow(sigma, 6.0);
158 # set the parameters for the appropriate type
159 self.cpp_force.setParams(globals.system_definition.getParticleData().getTypeByName(particle_type), lj1, lj2);
161 # track which particle types we have set
162 if not particle_type in self.particle_types_set:
163 self.particle_types_set.append(particle_type);
165 def update_coeffs(self):
166 # get a list of all particle types in the simulation
167 ntypes = globals.system_definition.getParticleData().getNTypes();
168 type_list = [];
169 for i in range(0,ntypes):
170 type_list.append(globals.system_definition.getParticleData().getNameByType(i));
172 # check to see if all particle types have been set
173 for cur_type in type_list:
174 if not cur_type in self.particle_types_set:
175 globals.msg.error(str(cur_type) + " coefficients missing in wall.lj\n");
176 raise RuntimeError("Error updating coefficients");