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: akohlmey
53 #pragma warning( push )
54 #pragma warning( disable : 4103 4244 )
57 #include <boost/python.hpp>
58 using namespace boost::python
;
60 #include "CGCMMForceCompute.h"
63 /*! \file CGCMMForceCompute.cc
64 \brief Defines the CGCMMForceCompute class
69 /*! \param sysdef System to compute forces on
70 \param nlist Neighborlist to use for computing the forces
71 \param r_cut Cuttoff radius beyond which the force is 0
72 \post memory is allocated and all parameters ljX are set to 0.0
74 CGCMMForceCompute::CGCMMForceCompute(boost::shared_ptr
<SystemDefinition
> sysdef
,
75 boost::shared_ptr
<NeighborList
> nlist
,
77 : ForceCompute(sysdef
), m_nlist(nlist
), m_r_cut(r_cut
)
79 m_exec_conf
->msg
->notice(5) << "Constructing CGCMMForceCompute" << endl
;
86 m_exec_conf
->msg
->error() << "pair.cgcmm: Negative r_cut makes no sense" << endl
;
87 throw runtime_error("Error initializing CGCMMForceCompute");
90 // initialize the number of types value
91 m_ntypes
= m_pdata
->getNTypes();
94 // allocate storage for lj12, lj9, lj6, and lj4 parameters
95 m_lj12
= new Scalar
[m_ntypes
*m_ntypes
];
96 m_lj9
= new Scalar
[m_ntypes
*m_ntypes
];
97 m_lj6
= new Scalar
[m_ntypes
*m_ntypes
];
98 m_lj4
= new Scalar
[m_ntypes
*m_ntypes
];
105 memset((void*)m_lj12
, 0, sizeof(Scalar
)*m_ntypes
*m_ntypes
);
106 memset((void*)m_lj9
, 0, sizeof(Scalar
)*m_ntypes
*m_ntypes
);
107 memset((void*)m_lj6
, 0, sizeof(Scalar
)*m_ntypes
*m_ntypes
);
108 memset((void*)m_lj4
, 0, sizeof(Scalar
)*m_ntypes
*m_ntypes
);
113 CGCMMForceCompute::~CGCMMForceCompute()
115 m_exec_conf
->msg
->notice(5) << "Destroying CGCMMForceCompute" << endl
;
117 // deallocate our memory
129 /*! \post The parameters \a lj12 through \a lj4 are set for the pairs \a typ1, \a typ2 and \a typ2, \a typ1.
130 \note \a lj? are low level parameters used in the calculation. In order to specify
131 these for a 12-4 and 9-6 lennard jones formula (with alpha), they should be set to the following.
134 - \a lj12 = 2.598076 * epsilon * pow(sigma,12.0)
137 - \a lj4 = -alpha * 2.598076 * epsilon * pow(sigma,4.0)
141 - \a lj9 = 6.75 * epsilon * pow(sigma,9.0);
142 - \a lj6 = -alpha * 6.75 * epsilon * pow(sigma,6.0)
146 - \a lj12 = 4.0 * epsilon * pow(sigma,12.0)
148 - \a lj6 = -alpha * 4.0 * epsilon * pow(sigma,4.0)
151 Setting the parameters for typ1,typ2 automatically sets the same parameters for typ2,typ1: there
152 is no need to call this funciton for symmetric pairs. Any pairs that this function is not called
153 for will have lj12 through lj4 set to 0.0.
155 \param typ1 Specifies one type of the pair
156 \param typ2 Specifies the second type of the pair
157 \param lj12 1/r^12 term
158 \param lj9 1/r^9 term
159 \param lj6 1/r^6 term
160 \param lj4 1/r^4 term
162 void CGCMMForceCompute::setParams(unsigned int typ1
, unsigned int typ2
, Scalar lj12
, Scalar lj9
, Scalar lj6
, Scalar lj4
)
164 if (typ1
>= m_ntypes
|| typ2
>= m_ntypes
)
166 m_exec_conf
->msg
->error() << "pair.cgcmm: Trying to set params for a non existant type! " << typ1
<< "," << typ2
<< endl
;
167 throw runtime_error("Error setting parameters in CGCMMForceCompute");
170 // set lj12 in both symmetric positions in the matrix
171 m_lj12
[typ1
*m_ntypes
+ typ2
] = lj12
;
172 m_lj12
[typ2
*m_ntypes
+ typ1
] = lj12
;
174 // set lj9 in both symmetric positions in the matrix
175 m_lj9
[typ1
*m_ntypes
+ typ2
] = lj9
;
176 m_lj9
[typ2
*m_ntypes
+ typ1
] = lj9
;
178 // set lj6 in both symmetric positions in the matrix
179 m_lj6
[typ1
*m_ntypes
+ typ2
] = lj6
;
180 m_lj6
[typ2
*m_ntypes
+ typ1
] = lj6
;
182 // set lj4 in both symmetric positions in the matrix
183 m_lj4
[typ1
*m_ntypes
+ typ2
] = lj4
;
184 m_lj4
[typ2
*m_ntypes
+ typ1
] = lj4
;
187 /*! CGCMMForceCompute provides
190 std::vector
< std::string
> CGCMMForceCompute::getProvidedLogQuantities()
193 list
.push_back("pair_cgcmm_energy");
197 Scalar
CGCMMForceCompute::getLogValue(const std::string
& quantity
, unsigned int timestep
)
199 if (quantity
== string("pair_cgcmm_energy"))
202 return calcEnergySum();
206 m_exec_conf
->msg
->error() << "pair.cgcmm: " << quantity
<< " is not a valid log quantity" << endl
;
207 throw runtime_error("Error getting log value");
211 /*! \post The CGCMM forces are computed for the given timestep. The neighborlist's
212 compute method is called to ensure that it is up to date.
214 \param timestep specifies the current time step of the simulation
216 void CGCMMForceCompute::computeForces(unsigned int timestep
)
218 // start by updating the neighborlist
219 m_nlist
->compute(timestep
);
221 // start the profile for this compute
222 if (m_prof
) m_prof
->push("CGCMM pair");
224 ArrayHandle
<Scalar4
> h_force(m_force
,access_location::host
, access_mode::overwrite
);
225 ArrayHandle
<Scalar
> h_virial(m_virial
,access_location::host
, access_mode::overwrite
);
226 const unsigned int virial_pitch
= m_virial
.getPitch();
228 // there are enough other checks on the input data: but it doesn't hurt to be safe
229 assert(h_force
.data
);
230 assert(h_virial
.data
);
232 // Zero data for force calculation.
233 memset((void*)h_force
.data
,0,sizeof(Scalar4
)*m_force
.getNumElements());
234 memset((void*)h_virial
.data
,0,sizeof(Scalar
)*m_virial
.getNumElements());
236 // depending on the neighborlist settings, we can take advantage of newton's third law
237 // to reduce computations at the cost of memory access complexity: set that flag now
238 bool third_law
= m_nlist
->getStorageMode() == NeighborList::half
;
240 // access the neighbor list
241 ArrayHandle
<unsigned int> h_n_neigh(m_nlist
->getNNeighArray(), access_location::host
, access_mode::read
);
242 ArrayHandle
<unsigned int> h_nlist(m_nlist
->getNListArray(), access_location::host
, access_mode::read
);
243 Index2D nli
= m_nlist
->getNListIndexer();
245 // access the particle data
246 ArrayHandle
< Scalar4
> h_pos(m_pdata
->getPositions(), access_location::host
, access_mode::read
);
248 assert(h_pos
.data
!= NULL
);
250 // get a local copy of the simulation box too
251 const BoxDim
& box
= m_pdata
->getBox();
253 // create a temporary copy of r_cut sqaured
254 Scalar r_cut_sq
= m_r_cut
* m_r_cut
;
256 // tally up the number of forces calculated
260 for (unsigned int i
= 0; i
< m_pdata
->getN(); i
++)
262 // access the particle's position and type (MEM TRANSFER: 4 scalars)
263 Scalar3 pi
= make_scalar3(h_pos
.data
[i
].x
, h_pos
.data
[i
].y
, h_pos
.data
[i
].z
);
264 unsigned int typei
= __scalar_as_int(h_pos
.data
[i
].w
);
266 assert(typei
< m_pdata
->getNTypes());
268 // access the lj12 and lj9 rows for the current particle type
269 Scalar
* __restrict__ lj12_row
= &(m_lj12
[typei
*m_ntypes
]);
270 Scalar
* __restrict__ lj9_row
= &(m_lj9
[typei
*m_ntypes
]);
271 Scalar
* __restrict__ lj6_row
= &(m_lj6
[typei
*m_ntypes
]);
272 Scalar
* __restrict__ lj4_row
= &(m_lj4
[typei
*m_ntypes
]);
274 // initialize current particle force, potential energy, and virial to 0
275 Scalar3 fi
= make_scalar3(0, 0, 0);
278 for (int k
= 0; k
< 6; k
++)
281 // loop over all of the neighbors of this particle
282 const unsigned int size
= (unsigned int)h_n_neigh
.data
[i
];
283 for (unsigned int j
= 0; j
< size
; j
++)
285 // increment our calculation counter
288 // access the index of this neighbor (MEM TRANSFER: 1 scalar)
289 unsigned int k
= h_nlist
.data
[nli(i
, j
)];
291 assert(k
< m_pdata
->getN());
293 // calculate dr (MEM TRANSFER: 3 scalars / FLOPS: 3)
294 Scalar3 pj
= make_scalar3(h_pos
.data
[k
].x
, h_pos
.data
[k
].y
, h_pos
.data
[k
].z
);
295 Scalar3 dx
= pi
- pj
;
297 // access the type of the neighbor particle (MEM TRANSFER: 1 scalar
298 unsigned int typej
= __scalar_as_int(h_pos
.data
[k
].w
);
300 assert(typej
< m_pdata
->getNTypes());
302 // apply periodic boundary conditions (FLOPS: 9 (worst case: first branch is missed, the 2nd is taken and the add is done)
303 dx
= box
.minImage(dx
);
305 // start computing the force
306 // calculate r squared (FLOPS: 5)
307 Scalar rsq
= dot(dx
, dx
);
309 // only compute the force if the particles are closer than the cuttoff (FLOPS: 1)
312 // compute the force magnitude/r in forcemag_divr (FLOPS: 14)
313 Scalar r2inv
= Scalar(1.0)/rsq
;
314 Scalar r3inv
= r2inv
/ sqrt(rsq
);
315 Scalar r6inv
= r3inv
* r3inv
;
316 Scalar forcemag_divr
= r6inv
* (r2inv
* (Scalar(12.0)*lj12_row
[typej
]*r6inv
+ Scalar(9.0)*r3inv
*lj9_row
[typej
]
317 + Scalar(6.0)*lj6_row
[typej
]) + Scalar(4.0)*lj4_row
[typej
]);
319 // compute the pair energy and virial (FLOPS: 6)
320 Scalar pair_virial
[6];
321 pair_virial
[0] = Scalar(0.5) * dx
.x
* dx
.x
* forcemag_divr
;
322 pair_virial
[1] = Scalar(0.5) * dx
.x
* dx
.y
* forcemag_divr
;
323 pair_virial
[2] = Scalar(0.5) * dx
.x
* dx
.z
* forcemag_divr
;
324 pair_virial
[3] = Scalar(0.5) * dx
.y
* dx
.y
* forcemag_divr
;
325 pair_virial
[4] = Scalar(0.5) * dx
.y
* dx
.z
* forcemag_divr
;
326 pair_virial
[5] = Scalar(0.5) * dx
.z
* dx
.z
* forcemag_divr
;
328 Scalar pair_eng
= Scalar(0.5) * (r6inv
* (lj12_row
[typej
] * r6inv
+ lj9_row
[typej
] * r3inv
+ lj6_row
[typej
]) + lj4_row
[typej
] * r2inv
* r2inv
);
330 // add the force, potential energy and virial to the particle i
332 fi
+= dx
*forcemag_divr
;
334 for (unsigned int l
= 0; l
< 6; l
++)
335 viriali
[l
] += pair_virial
[l
];
337 // add the force to particle j if we are using the third law (MEM TRANSFER: 10 scalars / FLOPS: 8)
340 h_force
.data
[k
].x
-= dx
.x
*forcemag_divr
;
341 h_force
.data
[k
].y
-= dx
.y
*forcemag_divr
;
342 h_force
.data
[k
].z
-= dx
.z
*forcemag_divr
;
343 h_force
.data
[k
].w
+= pair_eng
;
344 for (unsigned int l
= 0; l
< 6; l
++)
345 h_virial
.data
[l
*virial_pitch
+k
] += pair_virial
[l
];
351 // finally, increment the force, potential energy and virial for particle i
352 // (MEM TRANSFER: 10 scalars / FLOPS: 5)
353 h_force
.data
[i
].x
+= fi
.x
;
354 h_force
.data
[i
].y
+= fi
.y
;
355 h_force
.data
[i
].z
+= fi
.z
;
356 h_force
.data
[i
].w
+= pei
;
357 for (int l
= 0; l
< 6; l
++)
358 h_virial
.data
[l
*virial_pitch
+i
] += viriali
[l
];
361 int64_t flops
= m_pdata
->getN() * 5 + n_calc
* (3+5+9+1+14+6+8);
362 if (third_law
) flops
+= n_calc
* 8;
363 int64_t mem_transfer
= m_pdata
->getN() * (5+4+10)*sizeof(Scalar
) + n_calc
* (1+3+1)*sizeof(Scalar
);
364 if (third_law
) mem_transfer
+= n_calc
*10*sizeof(Scalar
);
365 if (m_prof
) m_prof
->pop(flops
, mem_transfer
);
368 void export_CGCMMForceCompute()
370 class_
<CGCMMForceCompute
, boost::shared_ptr
<CGCMMForceCompute
>, bases
<ForceCompute
>, boost::noncopyable
>
371 ("CGCMMForceCompute", init
< boost::shared_ptr
<SystemDefinition
>, boost::shared_ptr
<NeighborList
>, Scalar
>())
372 .def("setParams", &CGCMMForceCompute::setParams
)
377 #pragma warning( pop )