Thanks to visit codestin.com
Credit goes to github.com

Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
79 changes: 79 additions & 0 deletions hoomd/md/BoxDeformer.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,79 @@
// Copyright (c) 2009-2025 The Regents of the University of Michigan.
// Part of HOOMD-blue, released under the BSD 3-Clause License.

/*! \file BoxDeformer.cc
\brief Declares the BoxDeformer class
*/

#include "BoxDeformer.h"

using namespace std;

namespace hoomd
{
namespace md
{
/** @param sysdef System definition containing the particle data this method acts on
@param deltaT Time step
*/
BoxDeformer::BoxDeformer(std::shared_ptr<SystemDefinition> sysdef, Scalar deltaT)
: m_sysdef(sysdef), m_pdata(sysdef->getParticleData()), m_deltaT(deltaT)
{
assert(m_pdata);
}

BoxDeformer::~BoxDeformer() { }

void BoxDeformer::setDeltaT(Scalar deltaT)
{
if (deltaT < 0.0)
throw std::domain_error("delta_t must be positive");
m_deltaT = deltaT;
}

void BoxDeformer::update(uint64_t timestep)
{
// Get the current box
BoxDim old_box = m_pdata->getGlobalBox();

// Compute new box (child class will determine the new box geometry)
BoxDim new_box = computeNewBox(timestep);

if (new_box != old_box)
{
// Set the new global box
m_pdata->setGlobalBox(new_box);

// Apply any post-deformation processing
postDeformationProcessing(old_box, new_box);
}
}

// Post deformation particle processing: PBC wrapping by default but child classes can add up
void BoxDeformer::postDeformationProcessing(const BoxDim& old_box, BoxDim& new_box)
{
ArrayHandle<Scalar4> h_pos(m_pdata->getPositions(),
access_location::host,
access_mode::readwrite);
ArrayHandle<Scalar4> h_vel(m_pdata->getVelocities(),
access_location::host,
access_mode::readwrite);
ArrayHandle<int3> h_image(m_pdata->getImages(), access_location::host, access_mode::readwrite);

for (unsigned int i = 0; i < m_pdata->getN(); i++)
{
new_box.wrap(h_pos.data[i], h_vel.data[i], h_image.data[i]);
}
}

namespace detail
{
void export_BoxDeformer(pybind11::module& m)
{
pybind11::class_<BoxDeformer, std::shared_ptr<BoxDeformer>>(m, "BoxDeformer")
.def(pybind11::init<std::shared_ptr<SystemDefinition>, Scalar>());
}

} // end namespace detail
} // end namespace md
} // end namespace hoomd
64 changes: 64 additions & 0 deletions hoomd/md/BoxDeformer.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,64 @@
// Copyright (c) 2009-2025 The Regents of the University of Michigan.
// Part of HOOMD-blue, released under the BSD 3-Clause License.

/*! \file BoxDeformer.h
\brief Declares a base class for all box deformers
*/

#ifdef __HIPCC__
#error This header cannot be compiled by nvcc
#endif

#ifndef __BOX_DEFORMER_H__
#define __BOX_DEFORMER_H__

#include "hoomd/BoxDim.h"
#include "hoomd/ParticleData.h"
#include "hoomd/SystemDefinition.h"

#include <memory>
#include <pybind11/pybind11.h>
#include <string>

namespace hoomd
{
namespace md
{
/** Perform box deformation.
Implement methods for simulation under box deformation.
*/
class PYBIND11_EXPORT BoxDeformer
{
public:
/// Constructor
BoxDeformer(std::shared_ptr<SystemDefinition> sysdef, Scalar deltaT);

/// Destructor
virtual ~BoxDeformer();

/// Set the time step size to be consistent with the Integrator's
virtual void setDeltaT(Scalar deltaT);

// Apply deformation
void update(uint64_t timestep);

protected:
std::shared_ptr<SystemDefinition>
m_sysdef; //!< The system definition this method is associated with
std::shared_ptr<ParticleData> m_pdata; //!< The particle data this method is associated with
Scalar m_deltaT; //!< The time step

virtual BoxDim computeNewBox(uint64_t timestep);

virtual void postDeformationProcessing(const BoxDim& old_box, BoxDim& new_box);
};

namespace detail
{
/// Export the BoxDeformer class to python
void export_BoxDeformer(pybind11::module& m);
} // end namespace detail
} // end namespace md
} // end namespace hoomd

#endif // #ifndef __BOX_DEFORMER_H__
4 changes: 4 additions & 0 deletions hoomd/md/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -13,6 +13,7 @@ set(_md_sources module-md.cc
AlchemyData.cc
BendingRigidityMeshForceCompute.cc
BondTablePotential.cc
BoxDeformer.cc
CommunicatorGrid.cc
ComputeThermo.cc
ComputeThermoHMA.cc
Expand All @@ -31,6 +32,7 @@ set(_md_sources module-md.cc
HelfrichMeshForceCompute.cc
IntegrationMethodTwoStep.cc
IntegratorTwoStep.cc
LeesEdwardsBoxDeformer.cc
ManifoldZCylinder.cc
ManifoldDiamond.cc
ManifoldEllipsoid.cc
Expand Down Expand Up @@ -82,6 +84,7 @@ set(_md_headers ActiveForceComputeGPU.h
BendingRigidityMeshForceCompute.h
BondTablePotentialGPU.h
BondTablePotential.h
BoxDeformer.h
CommunicatorGridGPU.h
CommunicatorGrid.h
ComputeThermoGPU.cuh
Expand Down Expand Up @@ -164,6 +167,7 @@ set(_md_headers ActiveForceComputeGPU.h
HelfrichMeshForceCompute.h
IntegrationMethodTwoStep.h
IntegratorTwoStep.h
LeesEdwardsBoxDeformer.h
ManifoldZCylinder.h
ManifoldDiamond.h
ManifoldEllipsoid.h
Expand Down
11 changes: 11 additions & 0 deletions hoomd/md/IntegratorTwoStep.cc
Original file line number Diff line number Diff line change
Expand Up @@ -55,6 +55,12 @@ void IntegratorTwoStep::update(uint64_t timestep)
// ensure that prepRun() has been called
assert(m_prepared);

// pass timestep to deformer
if (m_deformer)
{
m_deformer->update(timestep);
}

// perform the first step of the integration on all groups
for (auto& method : m_methods)
{
Expand Down Expand Up @@ -142,6 +148,10 @@ void IntegratorTwoStep::setDeltaT(Scalar deltaT)
{
m_rigid_bodies->setDeltaT(deltaT);
}
if (m_deformer)
{
m_deformer->setDeltaT(deltaT);
}
}

/*! \param group Group over which to count degrees of freedom.
Expand Down Expand Up @@ -404,6 +414,7 @@ void export_IntegratorTwoStep(pybind11::module& m)
.def(pybind11::init<std::shared_ptr<SystemDefinition>, Scalar>())
.def_property_readonly("methods", &IntegratorTwoStep::getIntegrationMethods)
.def_property("rigid", &IntegratorTwoStep::getRigid, &IntegratorTwoStep::setRigid)
.def_property("deformer", &IntegratorTwoStep::getDeformer, &IntegratorTwoStep::setDeformer)
.def_property("integrate_rotational_dof",
&IntegratorTwoStep::getIntegrateRotationalDOF,
&IntegratorTwoStep::setIntegrateRotationalDOF)
Expand Down
13 changes: 13 additions & 0 deletions hoomd/md/IntegratorTwoStep.h
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,7 @@
#include "IntegrationMethodTwoStep.h"
#include "hoomd/Integrator.h"

#include "BoxDeformer.h"
#include "ForceComposite.h"

#pragma once
Expand Down Expand Up @@ -107,6 +108,16 @@ class PYBIND11_EXPORT IntegratorTwoStep : public Integrator
m_rigid_bodies = new_rigid;
}

/// Getter and setter for accessing box deformers
std::shared_ptr<BoxDeformer> getDeformer()
{
return m_deformer;
}
virtual void setDeformer(std::shared_ptr<BoxDeformer> deformer)
{
m_deformer = deformer;
}

/// Validate method groups.
void validateGroups();

Expand All @@ -116,6 +127,8 @@ class PYBIND11_EXPORT IntegratorTwoStep : public Integrator

std::shared_ptr<ForceComposite> m_rigid_bodies; /// definition and updater for rigid bodies

std::shared_ptr<BoxDeformer> m_deformer; /// box deformation methods

bool m_prepared; //!< True if preprun has been called

/// True when orientation degrees of freedom should be integrated
Expand Down
83 changes: 83 additions & 0 deletions hoomd/md/LeesEdwardsBoxDeformer.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,83 @@
// Copyright (c) 2009-2025 The Regents of the University of Michigan.
// Part of HOOMD-blue, released under the BSD 3-Clause License.

/*! \file LeesEdwardsBoxDeformer.cc
\brief Lees–Edwards box deformer for triclinic boxes.
*/

#include "LeesEdwardsBoxDeformer.h"

namespace hoomd
{
namespace md
{
/** @param sysdef System definition containing the particle data this method acts on
@param deltaT Time step
*/
LeesEdwardsBoxDeformer::LeesEdwardsBoxDeformer(std::shared_ptr<SystemDefinition> sysdef,
Scalar deltaT)
: BoxDeformer(sysdef, deltaT), m_new_box(m_pdata->getGlobalBox()),
m_xy(m_new_box.getTiltFactorXY()), m_xz(m_new_box.getTiltFactorXZ()),
m_yz(m_new_box.getTiltFactorYZ()), m_xy_rate(m_new_box.getTiltDeformationRateXY())
{
}

LeesEdwardsBoxDeformer::~LeesEdwardsBoxDeformer() { }

BoxDim LeesEdwardsBoxDeformer::computeNewBox(uint64_t timestep)
{
// Update tilt factor using stored rate and deltaT
m_xy += m_xy_rate * m_deltaT;

// Return updated BoxDim
m_new_box.setTiltFactors(m_xy, m_xz, m_yz);
m_new_box.setTiltDeformationRates(m_xy_rate, 0.0, 0.0);

return m_new_box;
}

void LeesEdwardsBoxDeformer::postDeformationProcessing(const BoxDim& old_box, BoxDim& new_box)
{
// Call base class to perform default PBC wrapping
BoxDeformer::postDeformationProcessing(old_box, new_box);

// Extra processing: box flipping and particle remapping
int flip = static_cast<int>(std::floor(m_xy + Scalar(0.5)));

if (flip != 0)
{
// Update stored tilt
m_xy -= Scalar(flip);
new_box.setTiltFactors(m_xy, m_xz, m_yz);

// Remap particle positions
ArrayHandle<Scalar4> h_pos(m_pdata->getPositions(),
access_location::host,
access_mode::readwrite);

Scalar Ly = new_box.getL().y;
for (unsigned int i = 0; i < m_pdata->getN(); i++)
{
h_pos.data[i].x -= Scalar(flip) * Ly;
}

m_new_box = new_box;
}
}

namespace detail
{
void export_LeesEdwardsBoxDeformer(pybind11::module& m)
{
pybind11::class_<LeesEdwardsBoxDeformer, BoxDeformer, std::shared_ptr<LeesEdwardsBoxDeformer>>(
m,
"LeesEdwardsBoxDeformer")
.def(pybind11::init<std::shared_ptr<SystemDefinition>, Scalar>())
.def_property("shear_rate",
&LeesEdwardsBoxDeformer::getShearRate,
&LeesEdwardsBoxDeformer::setShearRate);
}

} // end namespace detail
} // end namespace md
} // end namespace hoomd
66 changes: 66 additions & 0 deletions hoomd/md/LeesEdwardsBoxDeformer.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,66 @@
// Copyright (c) 2009-2025 The Regents of the University of Michigan.
// Part of HOOMD-blue, released under the BSD 3-Clause License.

/*! \file LeesEdwardsBoxDeformer.h
\brief Lees–Edwards box deformer for triclinic boxes.
*/

#ifdef __HIPCC__
#error This header cannot be compiled by nvcc
#endif

#ifndef __LEES_EDWARDS_BOX_DEFORMER_H__
#define __LEES_EDWARDS_BOX_DEFORMER_H__

#include "BoxDeformer.h"

namespace hoomd
{
namespace md
{
/*
Lees-Edwards deformation. Box flipping and particle remapping performed
when tilt exceeds 0.5.
*/
class PYBIND11_EXPORT LeesEdwardsBoxDeformer : public BoxDeformer
{
public:
LeesEdwardsBoxDeformer(std::shared_ptr<SystemDefinition> sysdef, Scalar deltaT);

virtual ~LeesEdwardsBoxDeformer();

/// Set the shear rate (in xy)
virtual void setShearRate(const Scalar xy_rate)
{
m_xy_rate = xy_rate;
}

/// Get the shear rate (in xy)
Scalar getShearRate()
{
return m_xy_rate;
}

protected:
BoxDim m_new_box; //!< box object
Scalar m_xy; //!< xy tilt
const Scalar m_xz; //!< xz tilt
const Scalar m_yz; //!< yz tilt
Scalar m_xy_rate; //!< shear rate, d(xy)/dt

/// Compute the new box based on the shear rate
BoxDim computeNewBox(uint64_t timestep) override;

/// Box flip and particle remapping (called after default PBC wrapping)
void postDeformationProcessing(const BoxDim& old_box, BoxDim& new_box) override;
};

namespace detail
{
/// Export LeesEdwardsBoxDeformer to python
void export_LeesEdwardsBoxDeformer(pybind11::module& m);
} // end namespace detail
} // end namespace md
} // end namespace hoomd

#endif
Loading
Loading