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
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.
50 // Maintainer: phillicl
52 #include <boost/shared_ptr.hpp>
54 #include "ForceCompute.h"
55 #include "BondedGroupData.h"
59 /*! \file TableAngleForceCompute.h
60 \brief Declares the TableAngleForceCompute class
64 #error This header cannot be compiled by nvcc
67 #ifndef __TABLEANGLEFORCECOMPUTE_H__
68 #define __TABLEANGLEFORCECOMPUTE_H__
70 //! Computes the potential and force on bonds based on values given in a table
72 Bond potentials and forces are evaluated for all bonded particle pairs in the system.
73 Both the potentials and forces are provided the tables V(r) and F(r) at discreet \a r values between \a thmin and
74 \a thmax. Evaluations are performed by simple linear interpolation, thus why F(r) must be explicitly specified to
75 avoid large errors resulting from the numerical derivative. Note that F(r) should store - dV/dr.
77 \b Table memory layout
79 V(\theta) and T(\theta) are specified for each bond type.
81 Three parameters need to be stored for each bond potential: thmin, thmax, and dr, the minimum r, maximum r, and spacing
82 between r values in the table respectively. For simple access on the GPU, these will be stored in a Scalar4 where
83 x is thmin, y is thmax, and z is dr.
85 V(0) is the value of V at r=thmin. V(i) is the value of V at r=thmin + dr * i where i is chosen such that r >= thmin
86 and r <= thmax. V(r) for r < thmin and > thmax is 0. The same goes for F. Thus V and F are defined between the region
87 [thmin,thmax], inclusive.
89 For ease of storing the data, all tables must be of the same number of points for all bonds.
92 Values are interpolated linearly between two points straddling the given r. For a given r, the first point needed, i
93 can be calculated via i = floorf((r - thmin) / dr). The fraction between ri and ri+1 can be calculated via
94 f = (r - thmin) / dr - Scalar(i). And the linear interpolation can then be performed via V(r) ~= Vi + f * (Vi+1 - Vi)
97 class TableAngleForceCompute
: public ForceCompute
100 //! Constructs the compute
101 TableAngleForceCompute(boost::shared_ptr
<SystemDefinition
> sysdef
,
102 unsigned int table_width
,
103 const std::string
& log_suffix
="");
106 virtual ~TableAngleForceCompute();
108 //! Set the table for a given type pair
109 virtual void setTable(unsigned int type
,
110 const std::vector
<Scalar
> &V
,
111 const std::vector
<Scalar
> &T
114 //! Returns a list of log quantities this compute calculates
115 virtual std::vector
< std::string
> getProvidedLogQuantities();
117 //! Calculates the requested log value and returns it
118 virtual Scalar
getLogValue(const std::string
& quantity
, unsigned int timestep
);
121 //! Get ghost particle fields requested by this pair potential
122 /*! \param timestep Current time step
124 virtual CommFlags
getRequestedCommFlags(unsigned int timestep
)
126 CommFlags flags
= CommFlags(0);
127 flags
[comm_flag::tag
] = 1;
128 flags
|= ForceCompute::getRequestedCommFlags(timestep
);
135 boost::shared_ptr
<AngleData
> m_angle_data
; //!< Angle data to use in computing angles
136 unsigned int m_table_width
; //!< Width of the tables in memory
137 GPUArray
<Scalar2
> m_tables
; //!< Stored V and T tables
138 Index2D m_table_value
; //!< Index table helper
139 std::string m_log_name
; //!< Cached log name
141 //! Actually compute the forces
142 virtual void computeForces(unsigned int timestep
);
145 //! Exports the TableAngleForceCompute class to python
146 void export_TableAngleForceCompute();