OpenMPCD
GaussianDumbbells.hpp
Go to the documentation of this file.
1 /**
2  * @file
3  * Defines the OpenMPCD::CUDA::MPCFluid::GaussianDumbbells class.
4  */
5 
6 #ifndef OPENMPCD_CUDA_MPCFLUID_GAUSSIANDUMBBELLS_HPP
7 #define OPENMPCD_CUDA_MPCFLUID_GAUSSIANDUMBBELLS_HPP
8 
10 
11 namespace OpenMPCD
12 {
13 namespace CUDA
14 {
15 namespace MPCFluid
16 {
17 
18  /**
19  * Fluid consisting of Gaussian dumbbells.
20  *
21  * The two particles interact with the following potential:
22  * \f[
23  * V
24  * = \frac{K}{2} \left( \vec{r}_1 - \vec{r}_2 \right)^2
25  * \f]
26  * where \f$ \vec{r}_i \f$ is the current position of particle \f$ i \f$
27  * and \f$ K \f$ is the spring constant.
28  *
29  * The configuration group of this fluid is expected to be named
30  * `mpc.fluid.dumbbell`, and contains:
31  * - The boolean value `analyticalStreaming`, which controls whether to
32  * use the analytically known solution to the equations of motion (value
33  * `true`), or whether to integrate the equations of motion using
34  * molecular dynamics (MD) during the streaming step (value `false`);
35  * - The floating-point value `rootMeanSquareLength`, which defines the
36  * root of the average of the squared bond length;
37  * - The floating-point value `zeroShearRelaxationTime`, which defines
38  * the relaxation time at zero shear.
39  * - If `analyticalStreaming` is `false`, the positive integer value
40  * `mdStepCount`, which specifies how many MD steps should be performed
41  * per MPC streaming step.
42  *
43  * The Gaussian dumbbell model is the one described in
44  * "Multiparticle collision dynamics simulations of viscoelastic fluids: Shear-thinning
45  * Gaussian dumbbells" by Bartosz Kowalik and Roland G. Winkler.
46  * Journal of Chemical Physics 138, 104903 (2013). http://dx.doi.org/10.1063/1.4792196
47  */
48  class GaussianDumbbells : public Base
49  {
50  public:
51  /**
52  * The constructor.
53  * @param[in] sim The simulation instance.
54  * @param[in] count The number of fluid particles.
55  * @param[in] streamingTimestep_ The timestep for a streaming step.
56  * @param[in] rng_ A random number generator to seed this instance's RNG with.
57  * @param[in] devMemMgr The Device memory manager.
58  * @param[in] kT The fluid temperature, times Boltzmann's constant.
59  * @param[in] leesEdwardsShearRate The shear rate for the Lees-Edwards boundary conditions.
60  */
61  GaussianDumbbells(const CUDA::Simulation* const sim, const unsigned int count,
62  const FP streamingTimestep_, RNG& rng_,
63  DeviceMemoryManager* const devMemMgr,
64  const FP kT, const FP leesEdwardsShearRate);
65 
66  /**
67  * The destructor.
68  */
70  {
71  }
72 
73  public:
74  virtual unsigned int getNumberOfLogicalEntities() const
75  {
76  return getParticleCount() / 2;
77  }
78 
80  {
81  return true;
82  }
83 
84  virtual unsigned int getNumberOfParticlesPerLogicalEntity() const
85  {
86  return 2;
87  }
88 
89  virtual void stream();
90 
91  private:
92  /**
93  * Reads the configuration.
94  */
95  void readConfiguration();
96 
97  /**
98  * Initializes the particle positions and velocities on the host.
99  * @param[in] leesEdwardsShearRate The shear rate for the Lees-Edwards boundary conditions.
100  */
101  void initializeOnHost(const FP leesEdwardsShearRate);
102 
103  /**
104  * Returns the initial position of the given particle's dumbbell partner.
105  * @param[in] position1 The position of the first particle of the dumbbell.
106  * @param[in] Wi The Weissenberg number.
107  */
109  getInitialDumbbellPartnerPosition(
111  const FP Wi) const;
112 
113  private:
114  bool streamAnalyticallyFlag; ///< Whether to use the analytical equations of motion or MD simulation.
115 
116  FP dumbbellRootMeanSquareLength; ///< The root-mean-square of the dumbbell bond length.
117  FP zeroShearRelaxationTime; ///< The zero-shear relaxation time of the dumbbells.
118  FP reducedSpringConstant; /**< The spring constant for the spring connecting the dumbbell constituents,
119  divided by the mass of each one constituent.*/
120 
121  unsigned int mdStepCount; ///< The number of velocity-Verlet steps in each streaming step.
122  };
123 
124 } //namespace MPCFluid
125 } //namespace CUDA
126 } //namespace OpenMPCD
127 
128 #endif
OpenMPCD::CUDA::MPCFluid::Base
Base class for MPC fluids.
Definition: CUDA/MPCFluid/Base.hpp:40
OpenMPCD::RemotelyStoredVector
Represents a vector whose data is stored elsewhere.
Definition: RemotelyStoredVector.hpp:26
OpenMPCD::CUDA::MPCFluid::GaussianDumbbells::numberOfParticlesPerLogicalEntityIsConstant
virtual bool numberOfParticlesPerLogicalEntityIsConstant() const
Returns whether all logical entities consist of the same number of MPC particles.
Definition: GaussianDumbbells.hpp:79
OpenMPCD::CUDA::DeviceMemoryManager
Class for managing memory on the CUDA Device.
Definition: DeviceMemoryManager.hpp:21
OpenMPCD::Vector3D
3-dimensional vector.
Definition: Vector3D.hpp:38
Base.hpp
OpenMPCD::CUDA::MPCFluid::GaussianDumbbells::GaussianDumbbells
GaussianDumbbells(const CUDA::Simulation *const sim, const unsigned int count, const FP streamingTimestep_, RNG &rng_, DeviceMemoryManager *const devMemMgr, const FP kT, const FP leesEdwardsShearRate)
The constructor.
OpenMPCD::CUDA::MPCFluid::GaussianDumbbells::~GaussianDumbbells
virtual ~GaussianDumbbells()
The destructor.
Definition: GaussianDumbbells.hpp:69
OpenMPCD::CUDA::MPCFluid::GaussianDumbbells::getNumberOfLogicalEntities
virtual unsigned int getNumberOfLogicalEntities() const
Returns the number of logical entities in the fluid.
Definition: GaussianDumbbells.hpp:74
OpenMPCD::CUDA::Simulation
MPCD simulation with Molecular Dynamics on CUDA-capable GPUs.
Definition: CUDA/Simulation.hpp:48
OpenMPCD::RNG
boost::mt11213b RNG
The random number generator type.
Definition: Types.hpp:18
OpenMPCD::CUDA::MPCFluid::GaussianDumbbells
Fluid consisting of Gaussian dumbbells.
Definition: GaussianDumbbells.hpp:48
OpenMPCD::CUDA::MPCFluid::Base::getParticleCount
unsigned int getParticleCount() const
Returns the number of MPC fluid particles.
Definition: CUDA/MPCFluid/Base.hpp:72
OpenMPCD::FP
double FP
Default floating point type.
Definition: Types.hpp:13
OpenMPCD::CUDA::MPCFluid::GaussianDumbbells::getNumberOfParticlesPerLogicalEntity
virtual unsigned int getNumberOfParticlesPerLogicalEntity() const
Returns the number of MPC particles per logical entity.
Definition: GaussianDumbbells.hpp:84
OpenMPCD::CUDA::MPCFluid::GaussianDumbbells::stream
virtual void stream()
Performs a streaming step.