OpenMPCD
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate > Class Template Referenceabstract

Class representing star polymers. More...

#include <StarPolymers.hpp>

Inheritance diagram for OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >:
Inheritance graph
[legend]
Collaboration diagram for OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >:
Collaboration graph
[legend]

Public Types

typedef PositionCoordinate ForceCoordinate
 The type to store force coordinates. More...
 
typedef ImplementationDetails::StarPolymers::ParticleType ParticleType
 Holds the enumeration of particle types. More...
 

Public Member Functions

 StarPolymers (const Configuration::Setting &settings, const BoundaryCondition::Base *const boundaryCondition)
 The constructor. More...
 
virtual ~StarPolymers ()
 The destructor. More...
 
void performMDTimestep ()
 Performs, on the Device, an MD timestep of size getMDTimeStepSize(). More...
 
std::size_t getStarCount () const
 Returns the total number of polymer stars in this instance. More...
 
std::size_t getArmCountPerStar () const
 Returns the number of arms per star. More...
 
std::size_t getParticleCountPerArm () const
 Returns the number of particles per arm. More...
 
bool hasMagneticParticles () const
 Returns whether arms have an additional magnetic particle attached. More...
 
std::size_t getParticleCountPerArmIncludingMagneticParticles () const
 Returns the number of particles per arm, including magnetic particles. More...
 
std::size_t getParticleCountPerStar () const
 Returns the total number of individual particles per star. More...
 
std::size_t getParticleCount () const
 Returns the total number of individual particles in all stars. More...
 
std::size_t getNumberOfLogicalEntities () const
 Returns the number of logical entities in the solute. More...
 
virtual FP getParticleMass () const
 Returns the mass of a particle, which is assumed to be equal for all particles in this instance. More...
 
void writeStructureToSnapshot (VTFSnapshotFile *const snapshot) const
 Writes structure information to the given snapshot file. More...
 
void writeStateToSnapshot (VTFSnapshotFile *const snapshot) const
 Writes the particle positions and velocities to the given snapshot file. More...
 
void readStateFromSnapshot (VTFSnapshotFile *const snapshot)
 Reads the next timestep block from the given snapshot file, and uses it as the current state. More...
 
bool snapshotStructureIsCompatible (const VTFSnapshotFile &snapshot) const
 Returns whether the given snapshot file is compatible with the star polymer architecture, e.g. More...
 
virtual std::size_t getParticleCount () const=0
 Returns the number of MPC solute particles. More...
 
virtual std::size_t getNumberOfLogicalEntities () const=0
 Returns the number of logical entities in the solute. More...
 
void fetchFromDevice () const
 Copies the MPC solute particles from the CUDA Device to the Host. More...
 
const RemotelyStoredVector< const MPCParticlePositionTypegetPosition (const std::size_t particleID) const
 Returns a MPC solute particle's position vector. More...
 
const RemotelyStoredVector< const MPCParticleVelocityTypegetVelocity (const std::size_t particleID) const
 Returns a MPC solute particle's velocity vector. More...
 
const MPCParticlePositionTypegetDevicePositions () const
 Returns a const pointer to the MPC solute positions on the Device. More...
 
MPCParticlePositionTypegetDevicePositions ()
 Returns a pointer to the MPC solute positions on the Device. More...
 
MPCParticlePositionTypegetHostPositions ()
 Returns a pointer to the MPC solute positions on the Host. More...
 
const MPCParticleVelocityTypegetDeviceVelocities () const
 Returns a const pointer to the MPC solute velocities on the Device. More...
 
MPCParticleVelocityTypegetDeviceVelocities ()
 Returns a pointer to the MPC solute velocities on the Device. More...
 
MPCParticleVelocityTypegetHostVelocities ()
 Returns a pointer to the MPC solute velocities on the Host. More...
 
FP getMDTimeStepSize () const
 Returns the number of time units advanced per MD step. More...
 
bool hasInstrumentation () const
 Returns whether the solute has instrumentation configured. More...
 
Instrumentation::BasegetInstrumentation () const
 Returns the solute instrumentation. More...
 
virtual FP getParticleMass () const=0
 Returns the mass of a particle, which is assumed to be equal for all particles in this instance. More...
 

Protected Member Functions

void pushToDevice ()
 Copies the MPC solute particles from the Host to the CUDA Device. More...
 

Protected Attributes

DeviceMemoryManager deviceMemoryManager
 The Device memory manager. More...
 
FP mdTimeStepSize
 The number of time units to advance per MD step. More...
 
Instrumentation::Baseinstrumentation
 The solute instrumentation. More...
 
boost::scoped_array< MPCParticlePositionTypeh_positions
 Host buffer for particle positions. More...
 
boost::scoped_array< MPCParticleVelocityTypeh_velocities
 Host buffer for particle velocities. More...
 
MPCParticlePositionTyped_positions
 Particle positions on the Device. More...
 
MPCParticleVelocityTyped_velocities
 Particle velocities on the Device. More...
 

Detailed Description

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
class OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >

Class representing star polymers.

A star polymer consists of one central core particle, attached to which there are a number of linear, homogeneous polymer chains, called "arms". The ends of the arms may be functionalized, by attaching an additional particle that carries a magnetic dipole moment.

The following three interaction potentials (magnetic, FENE, and WCA) are cumulative, i.e. a pair of particles can be subject to interaction due to two or more interactions at the same time.

The functionalized magnetic particles interact with one another via the OpenMPCD::PairPotentials::MagneticDipoleDipoleInteraction_ConstantIdenticalDipoles interaction. The prefactor property of the interaction class is set according to the interactionParameters.magneticPrefactor configuration property. The orientation of the dipoles may be set via the interactionParameters.dipoleOrientation setting, which must be an array of three floating-point values which define an axis. The vector formed by these three values need not be normalized; its magnitude is of no relevance.

Bonded particles – i.e. next neighbors in an arm, the central core particle and the particles attached to it, or the functionalized particles with their direct neighbor in the corresponding arm – interact via the OpenMPCD::PairPotentials::FENE potential. The interaction parameters K, l_0 and R are set to \( K_{\mu\nu} \), \( l_{\mu\nu} \), and \( R_{\mu\nu} \), respectively, which will be described below. The indices \( \mu \) and \( \nu \) denote between which kinds of particles (i.e. core, arm, or magnetic) the interaction effects a force.

Finally, each pair of particles interacts via the OpenMPCD::PairPotentials::WeeksChandlerAndersen_DistanceOffset (WCA) potential. The interaction parameters epsilon, sigma, and d are set to \( \varepsilon_{\mu\nu} \), \( \sigma_{\mu\nu} \), and \( D_{\mu\nu} \), respectively, which will be described below. The indices \( \mu \) and \( \nu \) have the same meaning as with the FENE potential.

In the configuration group interactionParameters, one has to set epsilon_core, epsilon_arm, and epsilon_magnetic, which will be referred to as \( \varepsilon_C \), \( \varepsilon_A \), and \( \varepsilon_M \), respectively. Similarly, the properties sigma_core, sigma_arm, and sigma_magnetic define \( \sigma_C \), \( \sigma_A \), and \( \sigma_M \), and D_core, D_arm, and D_magnetic define \( D_C \), \( D_A \), and \( D_M \). Finally, one also has to set magneticPrefactor in the interactionParameters group. The dipoleOrientation array, if not set, will be assumed to be [0, 0, 1].

From these quantities, the values for \( \varepsilon_{\mu\nu} \), \( \sigma{\mu\nu} \), and \( D_{\mu\nu} \), where each \( \mu \) and \( \nu \) can take on the values \( C \) for core particles, \( A \) for non-functionalized arm particles, and \( M \) for functionalized magnetic particles, can be calculated according to the Lorentz-Berthelot mixing rules:

\[ \varepsilon_{\mu\nu} = \sqrt{ \varepsilon_\mu \varepsilon_\nu } \]

\[ \sigma_{\mu\nu} = \frac{ \sigma_\mu + \sigma_\nu }{ 2 } \]

\[ D_{\mu\nu} = \frac{ D_\mu + D_\nu }{ 2 } \]

Finally, the parameters \( l_{\mu\nu} \), \( R_{\mu\nu} \), and \( K_{\mu\nu} \), are defined via

\[ l_{\mu\nu} = D_{\mu\nu} \]

\[ R_{\mu\nu} = 1.5 \sigma_{\mu\nu} \]

\[ K_{\mu\nu} = 30 \varepsilon_{\mu\nu} \sigma_{\mu\nu}^{-2} \]

In addition to the interactionParmeters settings group, the configuration passed to the constructor is expected to contain the structure group, which contains in starCount the number of stars, in armCountPerStar how many arms each star has, in armParticlesPerArm how many ordinary arm particles there are in each arm, and in hasMagneticParticles whether each arm should be terminated with a magnetic particle (in addition to the ordinary arm particles). Furthermore, particleMass contains the mass each individual particle has.

Furthermore, one must set mdTimeStepSize to the number of MPC time units one should advance per MD time-step, which must be a positive number.

If the configuration contains the value initialState, its value is taken as a path of a VTF snapshot file (see OpenMPCD::VTFSnapshotFile) that specifies the initial configuration of the stars. If the snapshot does not contain velocity information, initial velocities are assumed to be zero.

Template Parameters
PositionCoordinateThe type to store position coordinates.
VelocityCoordinateThe type to store velocity coordinates.

Definition at line 123 of file StarPolymers.hpp.

Member Typedef Documentation

◆ ForceCoordinate

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
typedef PositionCoordinate OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::ForceCoordinate

The type to store force coordinates.

Definition at line 127 of file StarPolymers.hpp.

◆ ParticleType

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
typedef ImplementationDetails::StarPolymers::ParticleType OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::ParticleType

Holds the enumeration of particle types.

Definition at line 130 of file StarPolymers.hpp.

Constructor & Destructor Documentation

◆ StarPolymers()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::StarPolymers ( const Configuration::Setting settings,
const BoundaryCondition::Base *const  boundaryCondition 
)

The constructor.

Exceptions
OpenMPCD::InvalidConfigurationExceptionThrows if the configuration is not valid.
Parameters
[in]settingsThe configuration for this instance.
[in]boundaryConditionThe boundary conditions that have been configured. Currently, only LeesEdwards boundary conditions are supported.

◆ ~StarPolymers()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
virtual OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::~StarPolymers ( )
virtual

The destructor.

Member Function Documentation

◆ fetchFromDevice()

Copies the MPC solute particles from the CUDA Device to the Host.

◆ getArmCountPerStar()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
std::size_t OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::getArmCountPerStar ( ) const

Returns the number of arms per star.

◆ getDevicePositions() [1/2]

Returns a pointer to the MPC solute positions on the Device.

Definition at line 116 of file CUDA/MPCSolute/Base.hpp.

◆ getDevicePositions() [2/2]

Returns a const pointer to the MPC solute positions on the Device.

Definition at line 108 of file CUDA/MPCSolute/Base.hpp.

◆ getDeviceVelocities() [1/2]

Returns a pointer to the MPC solute velocities on the Device.

Definition at line 140 of file CUDA/MPCSolute/Base.hpp.

◆ getDeviceVelocities() [2/2]

Returns a const pointer to the MPC solute velocities on the Device.

Definition at line 132 of file CUDA/MPCSolute/Base.hpp.

◆ getHostPositions()

Returns a pointer to the MPC solute positions on the Host.

Definition at line 124 of file CUDA/MPCSolute/Base.hpp.

◆ getHostVelocities()

Returns a pointer to the MPC solute velocities on the Host.

Definition at line 148 of file CUDA/MPCSolute/Base.hpp.

◆ getInstrumentation()

Returns the solute instrumentation.

Exceptions
OpenMPCD::NULLPointerExceptionIf OPENMPCD_DEBUG is defined, throws if !hasInstrumentation().

Definition at line 175 of file CUDA/MPCSolute/Base.hpp.

◆ getMDTimeStepSize()

Returns the number of time units advanced per MD step.

Definition at line 156 of file CUDA/MPCSolute/Base.hpp.

◆ getNumberOfLogicalEntities() [1/2]

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
std::size_t OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::getNumberOfLogicalEntities ( ) const
inline

Returns the number of logical entities in the solute.

Definition at line 206 of file StarPolymers.hpp.

◆ getNumberOfLogicalEntities() [2/2]

virtual std::size_t OpenMPCD::CUDA::MPCSolute::Base< MPCParticlePositionType , MPCParticleVelocityType >::getNumberOfLogicalEntities
pure virtualinherited

Returns the number of logical entities in the solute.

◆ getParticleCount() [1/2]

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
std::size_t OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::getParticleCount ( ) const

Returns the total number of individual particles in all stars.

◆ getParticleCount() [2/2]

virtual std::size_t OpenMPCD::CUDA::MPCSolute::Base< MPCParticlePositionType , MPCParticleVelocityType >::getParticleCount
pure virtualinherited

Returns the number of MPC solute particles.

◆ getParticleCountPerArm()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
std::size_t OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::getParticleCountPerArm ( ) const

Returns the number of particles per arm.

◆ getParticleCountPerArmIncludingMagneticParticles()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
std::size_t OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::getParticleCountPerArmIncludingMagneticParticles ( ) const

Returns the number of particles per arm, including magnetic particles.

◆ getParticleCountPerStar()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
std::size_t OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::getParticleCountPerStar ( ) const

Returns the total number of individual particles per star.

◆ getParticleMass() [1/2]

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
virtual FP OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::getParticleMass ( ) const
inlinevirtual

Returns the mass of a particle, which is assumed to be equal for all particles in this instance.

Definition at line 215 of file StarPolymers.hpp.

◆ getParticleMass() [2/2]

virtual FP OpenMPCD::CUDA::MPCSolute::Base< MPCParticlePositionType , MPCParticleVelocityType >::getParticleMass
pure virtualinherited

Returns the mass of a particle, which is assumed to be equal for all particles in this instance.

◆ getPosition()

const RemotelyStoredVector<const MPCParticlePositionType > OpenMPCD::CUDA::MPCSolute::Base< MPCParticlePositionType , MPCParticleVelocityType >::getPosition ( const std::size_t  particleID) const
inherited

Returns a MPC solute particle's position vector.

Warning
This function only returns the position that was current the last time fetchFromDevice was called.
Exceptions
OpenMPCD::OutOfBoundsExceptionIf OPENMPCD_DEBUG is defined, throws if particleID >= getParticleCount().
Parameters
[in]particleIDThe particle ID.

◆ getStarCount()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
std::size_t OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::getStarCount ( ) const
inline

Returns the total number of polymer stars in this instance.

Definition at line 168 of file StarPolymers.hpp.

◆ getVelocity()

const RemotelyStoredVector<const MPCParticleVelocityType > OpenMPCD::CUDA::MPCSolute::Base< MPCParticlePositionType , MPCParticleVelocityType >::getVelocity ( const std::size_t  particleID) const
inherited

Returns a MPC solute particle's velocity vector.

Warning
This function only returns the velocity that was current the last time fetchFromDevice was called.
Exceptions
OpenMPCD::OutOfBoundsExceptionIf OPENMPCD_DEBUG is defined, throws if particleID >= getParticleCount().
Parameters
[in]particleIDThe particle ID.

◆ hasInstrumentation()

Returns whether the solute has instrumentation configured.

Definition at line 164 of file CUDA/MPCSolute/Base.hpp.

◆ hasMagneticParticles()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
bool OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::hasMagneticParticles ( ) const

Returns whether arms have an additional magnetic particle attached.

◆ performMDTimestep()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
void OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::performMDTimestep ( )
virtual

◆ pushToDevice()

Copies the MPC solute particles from the Host to the CUDA Device.

◆ readStateFromSnapshot()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
void OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::readStateFromSnapshot ( VTFSnapshotFile *const  snapshot)

Reads the next timestep block from the given snapshot file, and uses it as the current state.

If no velocities are specified, the current ones are kept.

Exceptions
OpenMPCD::NULLPointerExceptionIf OPENMPCD_DEBUG is defined, throws if snapshot == nullptr.
OpenMPCD::InvalidArgumentExceptionThrows if !snapshotStructureIsCompatible().
Parameters
[in,out]snapshotThe snapshot file to read from.

◆ snapshotStructureIsCompatible()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
bool OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::snapshotStructureIsCompatible ( const VTFSnapshotFile snapshot) const

Returns whether the given snapshot file is compatible with the star polymer architecture, e.g.

whether the number of stars, arms per star, etc. match.

If !snapshot.isInReadMode(), false is returned.

Parameters
[in]snapshotThe snapshot to check for compatibility.

◆ writeStateToSnapshot()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
void OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::writeStateToSnapshot ( VTFSnapshotFile *const  snapshot) const

Writes the particle positions and velocities to the given snapshot file.

This function will call fetchFromDevice, and write the current Device data to the given snapshot file.

Exceptions
OpenMPCD::NULLPointerExceptionThrows if snapshot == nullptr.
OpenMPCD::InvalidArgumentExceptionThrows if the number of atoms declared in the snapshot does not match the number of particles in this instance.
OpenMPCD::InvalidArgumentExceptionThrows if the given snapshot is not in write mode.
Parameters
[in,out]snapshotThe snapshot file.

◆ writeStructureToSnapshot()

template<typename PositionCoordinate = MPCParticlePositionType, typename VelocityCoordinate = MPCParticleVelocityType>
void OpenMPCD::CUDA::MPCSolute::StarPolymers< PositionCoordinate, VelocityCoordinate >::writeStructureToSnapshot ( VTFSnapshotFile *const  snapshot) const

Writes structure information to the given snapshot file.

Exceptions
OpenMPCD::NULLPointerExceptionThrows if snapshot == nullptr.
OpenMPCD::InvalidArgumentExceptionThrows if the given snapshot is not in write mode, or if the structure block already has been processed.
Parameters
[in,out]snapshotThe snapshot file.

Member Data Documentation

◆ d_positions

Particle positions on the Device.

Definition at line 210 of file CUDA/MPCSolute/Base.hpp.

◆ d_velocities

Particle velocities on the Device.

Definition at line 211 of file CUDA/MPCSolute/Base.hpp.

◆ deviceMemoryManager

The Device memory manager.

Definition at line 199 of file CUDA/MPCSolute/Base.hpp.

◆ h_positions

boost::scoped_array<MPCParticlePositionType > OpenMPCD::CUDA::MPCSolute::Base< MPCParticlePositionType , MPCParticleVelocityType >::h_positions
mutableprotectedinherited

Host buffer for particle positions.

Definition at line 205 of file CUDA/MPCSolute/Base.hpp.

◆ h_velocities

boost::scoped_array<MPCParticleVelocityType > OpenMPCD::CUDA::MPCSolute::Base< MPCParticlePositionType , MPCParticleVelocityType >::h_velocities
mutableprotectedinherited

Host buffer for particle velocities.

Definition at line 207 of file CUDA/MPCSolute/Base.hpp.

◆ instrumentation

The solute instrumentation.

Definition at line 203 of file CUDA/MPCSolute/Base.hpp.

◆ mdTimeStepSize

The number of time units to advance per MD step.

Definition at line 201 of file CUDA/MPCSolute/Base.hpp.


The documentation for this class was generated from the following file: