Path Integral Quantum Monte Carlo
Public Member Functions
GrapheneLUTPotential Class Reference

Returns van der Waals' potential between a helium adatom and a graphene sheet using summation in reciprocal space. More...

#include <potential.h>

+ Inheritance diagram for GrapheneLUTPotential:
+ Collaboration diagram for GrapheneLUTPotential:

Public Member Functions

 GrapheneLUTPotential (const double, const double, const double, const double, const double, const Container *)
 Constructor.
 
 ~GrapheneLUTPotential ()
 Destructor.
 
double V (const dVec &r)
 Return the value of the van der Waals' interaction between a graphene sheet and a helium adatom at a position, r, above the sheet. More...
 
dVec gradV (const dVec &r)
 Return the gradient of the van der Waals' interaction between a graphene sheet and a helium adatom at a position, r, above the sheet. More...
 
Array< dVec, 1 > initialConfig (const Container *, MTRand &, const int)
 Initial configuration corresponding to graphene-helium vdW 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 grad2V (const dVec &)
 Grad^2 of the potential.
 
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...
 

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.
 

Detailed Description

Returns van der Waals' potential between a helium adatom and a graphene sheet using summation in reciprocal space.

Author: Nathan Nichols & Adrian Del Maestro Returns the potential energy resulting from a van der Waals' interaction between a helium adatom and a fixed infinite graphene lattice

Definition at line 1222 of file potential.h.

Member Function Documentation

◆ gradV()

dVec GrapheneLUTPotential::gradV ( const dVec r)
virtual

Return the gradient of the van der Waals' interaction between a graphene sheet and a helium adatom at a position, r, above the sheet.

Parameters
rthe position of a helium particle
Returns
the gradient of the van der Waals' potential for graphene-helium

Reimplemented from PotentialBase.

Definition at line 2902 of file potential.cpp.

2902  {
2903 
2904  return 0.0;
2905 
2906  /* double x = r[0]; */
2907  /* double y = r[1]; */
2908  /* double z = r[NDIM-1]+(Lzo2); */
2909  /* int zindex = int((z-zmin)/dr); */
2910  /* dVec dv = 0.0; */
2911 
2912  /* if (z < zmin) */
2913  /* return dv; */
2914 
2915  /* if (zindex >= tableLength) { */
2916  /* dv[NDIM-1] += q*(epsilon*sigma*sigma*2.*M_PI/A)* 4. * (pow(sigma/z,4) - pow(sigma/z,10)) / z; */
2917  /* return dv; */
2918  /* } */
2919 
2920  /* dv[NDIM-1] = gradvg(0,zindex); */
2921  /* double gdotb1 = 0.; */
2922  /* double gdotb2 = 0.; */
2923 
2924  /* int k = 0; */
2925 
2926  /* /1* NB: this has not been optimized!!! *1/ */
2927  /* for (int m = -gnum; m < gnum+1; m++) { */
2928  /* for (int n = -gnum; n < gnum+1; n++) { */
2929  /* k = karr(m+gnum,n+gnum); */
2930  /* if((m != 0) or (n != 0)) { */
2931  /* gdotb1 = ((m*g1x + n*g2x)*(b1x+x)) + ((m*g1y + n*g2y)*(b1y+y)); */
2932  /* gdotb2 = ((m*g1x + n*g2x)*(b2x+x)) + ((m*g1y + n*g2y)*(b2y+y)); */
2933 
2934  /* dv[0] += -(g1x*m+g2x*n) * (sin(gdotb1) + sin(gdotb2))*vg(k,zindex); */
2935  /* dv[1] += -(g1y*m+g2y*n) * (sin(gdotb1) + sin(gdotb2))*vg(k,zindex); */
2936  /* dv[NDIM-1] += (cos(gdotb1)+cos(gdotb2))*gradvg(k,zindex); */
2937  /* } */
2938  /* } */
2939  /* } */
2940  /* return dv; */
2941 }

◆ initialConfig()

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

Initial configuration corresponding to graphene-helium vdW potential.

Return an initial particle configuration.

We create particles at random locations above the graphene sheet.

Reimplemented from PotentialBase.

Definition at line 2948 of file potential.cpp.

2949  {
2950  /* The particle configuration */
2951  Array<dVec,1> initialPos(numParticles);
2952  initialPos = 0.0;
2953 
2954  int Nlayer = round(boxPtr->side[0]*boxPtr->side[1]/a1x/a1y/4);
2955  int numX = round(boxPtr->side[0]/a1x/2);
2956  int numY = round(boxPtr->side[1]/a1y/2);
2957  double gridX = boxPtr->side[0]/numX;
2958  double gridY = boxPtr->side[1]/numY;
2959  int gridSize = 0;
2960  if (3*Nlayer < numParticles){
2961  gridSize = numParticles - 3*Nlayer;
2962  }
2963  int numInLayers = numParticles - gridSize;
2964  int layerNum = 0;
2965  int iShifted = 0;
2966  dVec pos;
2967  for (int i = 0; i < numInLayers; i++){
2968  layerNum = i/Nlayer;
2969  iShifted = i - (layerNum*Nlayer);
2970  pos[0] = ((iShifted % numX) + 0.5)*gridX - (boxPtr->side[0]/2);
2971  pos[1] = ((iShifted / numX) + 0.5)*gridY - (boxPtr->side[1]/2);
2972  pos[NDIM-1] = ((layerNum+1)*3.0)-(boxPtr->side[NDIM-1]/2);
2973  boxPtr->putInside (pos);
2974  initialPos(i) = pos;
2975  }
2976  if (gridSize) {
2977  int initCubePart = ceil(pow(1.0*gridSize,1.0/(1.0*NDIM)));
2978  dVec initSide;
2979  for (int i = 0; i < NDIM - 1; i++) {
2980  initSide[i] = boxPtr->side[i]/initCubePart;
2981  }
2982 
2983  initSide[NDIM-1] = (boxPtr->side[NDIM-1] - 12.0)/initCubePart;
2984 
2985  for (int n = 0; n < gridSize; n++) {
2986  pos[0] = (((n % initCubePart) + 0.5) * initSide[0]) - (boxPtr->side[0]/2);
2987  pos[1] = (((n / initCubePart) + 0.5) * initSide[1]) - (boxPtr->side[1]/2);
2988  pos[NDIM-1] = (((n / initCubePart / initCubePart) + 0.5) * initSide[NDIM-1]) - (boxPtr->side[NDIM-1]/2) + 12.0;
2989  boxPtr->putInside (pos);
2990  initialPos(n + numInLayers) = pos;
2991  }
2992  }
2993  return initialPos;
2994 }
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
TinyVector< double, NDIM > dVec
A NDIM-vector of type double.
Definition: common.h:111
+ Here is the call graph for this function:

◆ V()

double GrapheneLUTPotential::V ( const dVec r)
virtual

Return the value of the van der Waals' interaction between a graphene sheet and a helium adatom at a position, r, above the sheet.

Parameters
rthe position of a helium particle
Returns
the van der Waals' potential for graphene-helium

Reimplemented from PotentialBase.

Definition at line 2853 of file potential.cpp.

2853  {
2854 
2855  double z = r[2]+(Lzo2);
2856  if (z < zmin)
2857  return V_zmin;
2858 
2859  /* A sigmoid to represent the hard-wall */
2860  double VWall = V_zmin/(1.0+exp(-invWallWidth*(z-zWall)));
2861 
2862  int zindex = int((z-zmin)/dr);
2863  if (zindex >= tableLength)
2864  return q*(epsilon*sigma*sigma*2.*M_PI/A)*( ((2./5.)*pow((sigma/z),10)) - pow((sigma/z),4) );
2865 
2866  double x1 = r[0]+b1x;
2867  double x2 = r[0]+b2x;
2868  double y1 = r[1]+b1y;
2869  double y2 = r[1]+b2y;
2870 
2871  double mx,my;
2872  double v = vg(0,zindex);
2873 
2874  int ig = 1;
2875  for (int m = -gnum; m <= gnum; m++) {
2876  mx = m*g1x;
2877  my = m*g1y;
2878  for (int n = 1; n <= gnum; n++) {
2879  v += 2.0*(cos((mx + n*g2x)*x1 + (my + n*g2y)*y1) +
2880  cos((mx + n*g2x)*x2 + (my + n*g2y)*y2)) * vg(gMagID(ig),zindex);
2881  ig++;
2882  }
2883  }
2884 
2885  for (int m=1; m <= gnum; m++){
2886  mx = m*g1x;
2887  my = m*g1y;
2888  v += 2.0*(cos(mx*x1 + my*y1) + cos(mx*x2 + my*y2)) * vg(gMagID(ig),zindex);
2889  ig++;
2890  }
2891 
2892  return v + VWall;
2893 }

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