Path Integral Quantum Monte Carlo
Public Member Functions
PlatedLJCylinderPotential Class Reference

Computes the value of the external wall potential for a plated cylindrical cavity. More...

#include <potential.h>

+ Inheritance diagram for PlatedLJCylinderPotential:
+ Collaboration diagram for PlatedLJCylinderPotential:

Public Member Functions

 PlatedLJCylinderPotential (const double, const double, const double, const double, const double)
 Constructor. More...
 
 ~PlatedLJCylinderPotential ()
 Destructor.
 
double V (const dVec &r)
 The integrated LJ Wall potential.
 
dVec gradV (const dVec &)
 Return the gradient of aziz potential for separation r using a lookup table.
 
double grad2V (const dVec &)
 Return the Laplacian of aziz potential for separation r using a lookup table.
 
Array< dVec, 1 > initialConfig (const Container *, MTRand &, const int)
 Initial configuration corresponding to the LJ cylinder potential. More...
 
- Public Member Functions inherited from PotentialBase
 PotentialBase ()
 Constructor.
 
virtual ~PotentialBase ()
 Destructor.
 
virtual double V (const dVec &, const dVec &)
 The effective potential for the pair product approximation.
 
virtual double dVdlambda (const dVec &, const dVec &)
 The derivative of the effective potential with respect to lambda and tau.
 
virtual double dVdtau (const dVec &, const dVec &)
 
void output (const double)
 A debug method that output's the potential to a supplied separation. More...
 
virtual Array< double, 1 > getExcLen ()
 Array to hold data elements. More...
 
- Public Member Functions inherited from TabulatedPotential
 TabulatedPotential ()
 Constructor.
 
virtual ~TabulatedPotential ()
 Destructor.
 

Additional Inherited Members

- Data Fields inherited from PotentialBase
double tailV
 Tail correction factor.
 
- Protected Member Functions inherited from PotentialBase
double deltaSeparation (double sep1, double sep2) const
 Return the minimum image difference for 1D separations.
 
- Protected Member Functions inherited from TabulatedPotential
void initLookupTable (const double, const double)
 Given a discretization factor and the system size, create and fill the lookup tables for the potential and its derivative.
 
virtual double newtonGregory (const Array< double, 1 > &, const TinyVector< double, 2 > &, const double)
 Use the Newton-Gregory forward difference method to do a 2-point lookup on the potential table. More...
 
virtual double direct (const Array< double, 1 > &, const TinyVector< double, 2 > &, const double)
 Use a direct lookup for the potential table. More...
 
- Protected Attributes inherited from TabulatedPotential
Array< double, 1 > lookupV
 A potential lookup table.
 
Array< double, 1 > lookupdVdr
 A lookup table for dVint/dr.
 
Array< double, 1 > lookupd2Vdr2
 A lookup table for d2Vint/dr2.
 
double dr
 The discretization for the lookup table.
 
int tableLength
 The number of elements in the lookup table.
 
TinyVector< double, 2 > extV
 Extremal value of V.
 
TinyVector< double, 2 > extdVdr
 Extremal value of dV/dr.
 
TinyVector< double, 2 > extd2Vdr2
 Extremal value of d2V/dr2.
 

Detailed Description

Computes the value of the external wall potential for a plated cylindrical cavity.

Definition at line 440 of file potential.h.

Constructor & Destructor Documentation

◆ PlatedLJCylinderPotential()

PlatedLJCylinderPotential::PlatedLJCylinderPotential ( const double  Ro_,
const double  Rw,
const double  sigmaPlated_,
const double  epsilonPlated_,
const double  densityPlated_ 
)

Constructor.

We create a nested Lennard-Jones cylinder, which uses a lookup table to hold the value of the integrated 6-12 potential for helium atoms interacting with the walls of a nanpore pre-plated with an adsorbed gas.

Parameters
radiusThe radius of the cylinder

Definition at line 800 of file potential.cpp.

801  :
802  PotentialBase(),
804 {
805  /* The outer radius of the tube */
806  Ro = Ro_;
807 
808  /* Plating substrates */
809  densityPlated = densityPlated_; //0.207; // atoms / angstrom^3
810  epsilonPlated = epsilonPlated_; //36.13613; // Kelvin
811  sigmaPlated = sigmaPlated_; //3.0225; // angstroms
812  Ri = Ro-Rw; //3.5;
813 
814  // Neon Values
815  /* densityPlated = 0.207; // atoms / angstrom^3 FIXME */
816  /* epsilonPlated = 20.161242; // Kelvin */
817  /* sigmaPlated = 2.711; // angstroms */
818  /* Ri = Ro-(1.52*2); //FIXME Fix hardcoded R2 */
819 
820  /* Substrate parameters for MCM-41 obtained by fitting
821  * @see https://nano.delmaestro.org/index.php?title=Effective_external_potential_(nanopores)
822  */
823  density = 1.000; // atoms / angstrom^3
824  epsilon = 1.59; // Kelvin
825  sigma = 3.44; // angstroms
826 
827  /* We choose a mesh consisting of 10^6 points, and create the lookup table */
828  dR = (1.0E-6)*Ri;
829  initLookupTable(dR,Ri);
830 
831  /* Find the minimun of the potential */
832  minV = 1.0E5;
833  for (int n = 0; n < tableLength; n++) {
834  if (lookupV(n) < minV)
835  minV = lookupV(n);
836  }
837 
838  /* The extremal values for the lookup table */
839  extV = valueV(0.0),valueV(Ri);
840  extdVdr = valuedVdr(0.0),valuedVdr(Ri);
841 }
PotentialBase()
Constructor.
Definition: potential.cpp:25
TabulatedPotential()
Constructor.
Definition: potential.cpp:149
int tableLength
The number of elements in the lookup table.
Definition: potential.h:91
TinyVector< double, 2 > extdVdr
Extremal value of dV/dr.
Definition: potential.h:94
TinyVector< double, 2 > extV
Extremal value of V.
Definition: potential.h:93
void initLookupTable(const double, const double)
Given a discretization factor and the system size, create and fill the lookup tables for the potentia...
Definition: potential.cpp:168
Array< double, 1 > lookupV
A potential lookup table.
Definition: potential.h:86
+ Here is the call graph for this function:

Member Function Documentation

◆ initialConfig()

Array< dVec, 1 > PlatedLJCylinderPotential::initialConfig ( const Container boxPtr,
MTRand &  random,
const int  numParticles 
)
virtual

Initial configuration corresponding to the LJ cylinder potential.

Return an initial particle configuration.

Return a set of initial positions inside the cylinder.

Reimplemented from PotentialBase.

Definition at line 948 of file potential.cpp.

949  {
950 
951  /* The particle configuration */
952  Array<dVec,1> initialPos(numParticles);
953  initialPos = 0.0;
954 
955  /* We shift the radius inward to account for the excluded volume from the
956  * hard wall. This represents the largest prism that can be put
957  * inside a cylinder. */
958  dVec lside;
959  lside[0] = lside[1] = sqrt(2.0)*(Ri-sigmaPlated);
960  lside[2] = boxPtr->side[NDIM-1];
961 
962  /* Get the linear size per particle */
963  double initSide = pow((1.0*numParticles/product(lside)),-1.0/(1.0*NDIM));
964 
965  /* We determine the number of initial grid boxes there are in
966  * in each dimension and compute their size */
967  int totNumGridBoxes = 1;
968  iVec numNNGrid;
969  dVec sizeNNGrid;
970 
971  for (int i = 0; i < NDIM; i++) {
972  numNNGrid[i] = static_cast<int>(ceil((lside[i] / initSide) - EPS));
973 
974  /* Make sure we have at least one grid box */
975  if (numNNGrid[i] < 1)
976  numNNGrid[i] = 1;
977 
978  /* Compute the actual size of the grid */
979  sizeNNGrid[i] = lside[i] / (1.0 * numNNGrid[i]);
980 
981  /* Determine the total number of grid boxes */
982  totNumGridBoxes *= numNNGrid[i];
983  }
984 
985  /* Now, we place the particles at the middle of each box */
986  PIMC_ASSERT(totNumGridBoxes>=numParticles);
987  dVec pos;
988  for (int n = 0; n < totNumGridBoxes; n++) {
989 
990  iVec gridIndex;
991  for (int i = 0; i < NDIM; i++) {
992  int scale = 1;
993  for (int j = i+1; j < NDIM; j++)
994  scale *= numNNGrid[j];
995  gridIndex[i] = (n/scale) % numNNGrid[i];
996  }
997 
998  for (int i = 0; i < NDIM; i++)
999  pos[i] = (gridIndex[i]+0.5)*sizeNNGrid[i] - 0.5*lside[i];
1000 
1001  boxPtr->putInside(pos);
1002 
1003  if (n < numParticles)
1004  initialPos(n) = pos;
1005  else
1006  break;
1007  }
1008 
1009  return initialPos;
1010 }
virtual void putInside(dVec &) const =0
Place a vector inside the simulation cell.
dVec side
The linear dimensions of the box.
Definition: container.h:31
#define NDIM
Number of spatial dimnsions.
Definition: common.h:71
#define EPS
A small number.
Definition: common.h:94
TinyVector< double, NDIM > dVec
A NDIM-vector of type double.
Definition: common.h:111
#define PIMC_ASSERT(X)
Rename assert method.
Definition: common.h:64
TinyVector< int, NDIM > iVec
A NDIM-vector of type integer.
Definition: common.h:114
+ Here is the call graph for this function:

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