Path Integral Quantum Monte Carlo
Public Member Functions | Data Fields | Protected Member Functions | Protected Attributes
ActionBase Class Reference

Holds a base class that all action classes will be derived from. More...

#include <action.h>

+ Inheritance diagram for ActionBase:
+ Collaboration diagram for ActionBase:

Public Member Functions

 ActionBase (const Path &, LookupTable &, PotentialBase *, PotentialBase *, WaveFunctionBase *, bool _local=true, string _name="Base", double _endFactor=1.0, int _period=1)
 Setup the path data members and canonical re-weighting factors.
 
virtual ~ActionBase ()
 Empty base constructor.
 
string getActionName ()
 Returns the action name.
 
double kineticAction ()
 The full kinetic Action
More...
 
double kineticAction (const beadLocator &)
 The kinetic Action at a single slice
More...
 
double kineticAction (const beadLocator &, int wlLength)
 The kinetic Action for wlLength slices. More...
 
virtual double potentialAction ()
 The effective potential inter-ACTION for various pass conditions.
 
virtual double potentialAction (const beadLocator &, const beadLocator &)
 Return the potential action for a path. More...
 
virtual double potentialAction (const beadLocator &)
 
virtual double barePotentialAction (const beadLocator &)
 
virtual double potentialActionCorrection (const beadLocator &)
 
virtual double potentialActionCorrection (const beadLocator &, const beadLocator &)
 
virtual double derivPotentialActionTau (int)
 
virtual double derivPotentialActionLambda (int)
 
virtual double secondderivPotentialActionTau (int)
 
virtual double derivPotentialActionTau (int, double)
 
virtual double derivPotentialActionLambda (int, double)
 
virtual dVec gradPotentialAction (int)
 
virtual double rDOTgradUterm1 (int)
 
virtual double rDOTgradUterm2 (int)
 
virtual double deltaDOTgradUterm1 (int)
 
virtual double deltaDOTgradUterm2 (int)
 
virtual double virKinCorr (int)
 
virtual TinyVector< double, 2 > potential (int)
 
virtual double potential (int, double)
 
void setShift (int _shift)
 The public method that sets the tau scaling factor.
 
int getShift ()
 Get the tau scaling factor.
 
double rho0 (const dVec &, const dVec &, int)
 The free-particle density matrix. More...
 
double rho0 (const beadLocator &, const beadLocator &, int)
 The free-particle density matrix for two beadLocators with imaginary time separation M.
 
double rho0 (const dVec &, const int)
 The free-particle density matrix for a given spatial and temporal separation. More...
 
double ensembleWeight (const int)
 The ensemble particle number weighting factor. More...
 

Data Fields

const bool local
 Is the action local in imaginary time?
 
const int period
 
PotentialBaseexternalPtr
 The external potential.
 
PotentialBaseinteractionPtr
 The interaction potential.
 
Array< int, 1 > sepHist
 A histogram of separations.
 
Array< int, 1 > cylSepHist
 A histogram of separations for a cylinder.
 

Protected Member Functions

double tau ()
 The local shifted value of tau.
 
void updateSepHist (const dVec &)
 Update the separation histogram. More...
 

Protected Attributes

string name
 The name of the action.
 
LookupTablelookup
 We need a non-constant reference for updates.
 
const Pathpath
 A reference to the paths.
 
WaveFunctionBasewaveFunctionPtr
 A pointer to a trial wave function object.
 
double endFactor
 Mutiplictive factor of the potential action on ends.
 
int shift
 The scaling factor for tau.
 
bool canonical
 Are we in the canonical ensemble?
 
int numBeads0
 The target number of beads.
 
bool window
 
int windowWidth
 
bool gaussianEnsemble
 
double gaussianEnsembleSD
 
beadLocator bead2
 
beadLocator bead3
 
dVec sep
 
dVec sep2
 
double dSep
 

Detailed Description

Holds a base class that all action classes will be derived from.

Implements the details of the action, including the potential and kinetic pieces. Two types of actions, local and non-local derive from this base class and the actually used actions then derive from these.

Definition at line 29 of file action.h.

Member Function Documentation

◆ ensembleWeight()

double ActionBase::ensembleWeight ( const int  deltaNumBeads)

The ensemble particle number weighting factor.

The ensemble weighting factor.

For the grand canonical ensemble we return 1, for the canonical ensemble we return exp(-[(nb'-nb_0)^2 - (nb-nb_0)^2]/dN^2) where nb = number of beads

Definition at line 122 of file action.cpp.

122  {
123 
124  // ensembleWieght returns 1.0 unless window or ensemble weight is used
125  if ( window || gaussianEnsemble ){
126  int numBeads = path.worm.getNumBeadsOn();
127  int numBeadsP = numBeads + deltaNumBeads;
128 
129  // Check if we are within bead window
130  if (window){
131  int beadWindowWidth = path.numTimeSlices*windowWidth;
132  if ( ( numBeadsP > (numBeads0+beadWindowWidth) )||
133  (numBeadsP < (numBeads0-beadWindowWidth)) )
134  return 0.0;
135  }
136  // If we are within window check gaussian weight
137  if (gaussianEnsemble){
138  double sigmaBeads = path.numTimeSlices*gaussianEnsembleSD;
139  double xp = 1.0*(numBeadsP-numBeads0)*(numBeadsP-numBeads0)/
140  (2.0*sigmaBeads*sigmaBeads);
141  double x = 1.0*(numBeads-numBeads0)*(numBeads-numBeads0)/
142  (2.0*sigmaBeads*sigmaBeads);
143  return exp(-xp+x);
144  }
145  else
146  return 1.0;
147  }
148  else
149  return 1.0; // Return 1.0 if no enembleWeight is used
150 
151 }
int numBeads0
The target number of beads.
Definition: action.h:126
const Path & path
A reference to the paths.
Definition: action.h:118
const int numTimeSlices
A local constant copy of the number of time slices.
Definition: path.h:37
Worm worm
Details on the worm.
Definition: path.h:44
int getNumBeadsOn() const
Return the number of active beads.
Definition: worm.h:88
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ kineticAction() [1/3]

double ActionBase::kineticAction ( )

The full kinetic Action

Return the total kinetic action.

The total kinetic action evaluated for all particles and all time slices.

Definition at line 205 of file action.cpp.

205  {
206 
207  double totK = 0.0;
208 
209  /* Calculate the kinetic energy. Even though there
210  * may be multiple mixing and swaps, it doesn't matter as we always
211  * just advance one time step at a time, as taken care of through the
212  * linking arrays. This has been checked! */
213  beadLocator beadIndex;
214  for (int slice = 0; slice < path.numTimeSlices; slice++) {
215  for (int ptcl = 0; ptcl < path.numBeadsAtSlice(slice); ptcl++) {
216  beadIndex = slice,ptcl;
217  totK += kineticAction(beadIndex);
218  }
219  }
220 
221  return totK;
222 }
double kineticAction()
The full kinetic Action
Definition: action.cpp:205
Array< int, 1 > numBeadsAtSlice
The number of active beads at a given time slice.
Definition: path.h:48
TinyVector< int, 2 > beadLocator
time-slice,bead-number world line index
Definition: common.h:117
+ Here is the caller graph for this function:

◆ kineticAction() [2/3]

double ActionBase::kineticAction ( const beadLocator beadIndex)

The kinetic Action at a single slice

Return the kinetic action for a single bead.

The kinetic action involves two contributions as each bead is linked to two other ones at the next and previous time slice.

Definition at line 159 of file action.cpp.

159  {
160 
161  double totK = 0.0;
162  /* We get the Kinetic action for a single slice, which connects
163  * to two other slices */
164  dVec vel;
165  vel = path.getVelocity(path.prev(beadIndex));
166  totK += dot(vel,vel);
167  vel = path.getVelocity(beadIndex);
168  totK += dot(vel,vel);
169 
170  return totK / (4.0*constants()->lambda()*tau());
171 }
double tau()
The local shifted value of tau.
Definition: action.h:133
double lambda() const
Get lambda = hbar^2/(2mk_B)
Definition: constants.h:46
dVec getVelocity(const beadLocator &) const
Return the velocity between two time slices of a given particle as a ndim-vector.
Definition: path.h:183
beadLocator & prev(int slice, int ptcl)
Move one link backward in imaginary time.
Definition: path.h:95
TinyVector< double, NDIM > dVec
A NDIM-vector of type double.
Definition: common.h:111
ConstantParameters * constants()
Global public access to the constants.
Definition: constants.h:201
+ Here is the call graph for this function:

◆ kineticAction() [3/3]

double ActionBase::kineticAction ( const beadLocator beadIndex,
int  wlLength 
)

The kinetic Action for wlLength slices.

Return the kinetic action for wlLength time slices starting at the bead given by beadIndex.


The total kinetic action for a string of connected beads all on the same worldline.

Definition at line 180 of file action.cpp.

180  {
181 
182  double totK = 0.0;
183 
184  /* Make a local copy of beadIndex */
185  beadLocator beadID;
186  beadID = beadIndex;
187 
188  /* Starting from the initial bead, move wlLength slices in imaginary time
189  * and accumulate the kinetic action */
190  dVec vel;
191  for (int n = 0; n < wlLength; n++) {
192  vel = path.getVelocity(beadID);
193  totK += dot(vel, vel);
194  beadID = path.next(beadID);
195  }
196 
197  return totK / (4.0*constants()->lambda()*tau());
198 }
beadLocator & next(int slice, int ptcl)
Move one link forward in imaginary time.
Definition: path.h:86
+ Here is the call graph for this function:

◆ potentialAction()

double ActionBase::potentialAction ( const beadLocator startBead,
const beadLocator endBead 
)
virtual

Return the potential action for a path.

Given a starting and ending bead, compute the potential action for the path.

NOTE: if startBead == endBead this returns potentialAction(startBead) I'm not sure if this is the desired behavior yet.

Parameters
startBeadthe starting bead
endBeadthe end bead on the path
Returns
the total potential action for a given path

Definition at line 237 of file action.cpp.

238  {
239 
240  double totU = 0.0;
241  beadLocator beadIndex;
242  beadIndex = startBead;
243  do {
244  totU += potentialAction(beadIndex);
245  beadIndex = path.next(beadIndex);
246  } while (!all(beadIndex==path.next(endBead)));
247  return totU;
248 }
virtual double potentialAction()
The effective potential inter-ACTION for various pass conditions.
Definition: action.h:48
+ Here is the call graph for this function:

◆ rho0() [1/2]

double ActionBase::rho0 ( const dVec r1,
const dVec r2,
int  M 
)

The free-particle density matrix.

The free-particle density matrix for two particles with imaginary time separation M.

Definition at line 83 of file action.cpp.

83  {
84  dVec vel;
85  vel = r2 - r1;
86  path.boxPtr->putInBC(vel);
87  double rho0Norm = pow(4.0*M_PI*constants()->lambda()*M*constants()->tau(),-0.5*NDIM);
88  return ( rho0Norm * exp(-dot(vel,vel) / (4.0*constants()->lambda()*M*constants()->tau()) ) );
89 }
void putInBC(dVec &r) const
Place a vector in boundary conditions.
Definition: container.h:50
const Container * boxPtr
A constant reference to the container class.
Definition: path.h:43
#define NDIM
Number of spatial dimnsions.
Definition: common.h:71
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ rho0() [2/2]

double ActionBase::rho0 ( const dVec vel,
const int  M 
)

The free-particle density matrix for a given spatial and temporal separation.

Parameters
velThe spatial separtion R - R'
MThe imaginary time separation
Returns
The free density matrix

Definition at line 110 of file action.cpp.

110  {
111  double rho0Norm = pow(4.0*M_PI*constants()->lambda()*M*constants()->tau(),-0.5*NDIM);
112  return ( rho0Norm * exp(-dot(vel,vel) / (4.0*constants()->lambda()*M*constants()->tau()) ) );
113 }
+ Here is the call graph for this function:

◆ updateSepHist()

void ActionBase::updateSepHist ( const dVec _sep)
inlineprotected

Update the separation histogram.

We multiply by the periodic array, as we will only measure separations along spatial dimensions with periodic boundary conditions.

Definition at line 71 of file action.cpp.

71  {
72  dVec psep;
73  psep = path.boxPtr->periodic*_sep;
74  int nR = int(sqrt(dot(psep,psep))/dSep);
75  if (nR < NPCFSEP)
76  ++sepHist(nR);
77 }
Array< int, 1 > sepHist
A histogram of separations.
Definition: action.h:111
TinyVector< unsigned int, NDIM > periodic
Determines which dimensions have periodic bc.
Definition: container.h:29
#define NPCFSEP
Spatial separations to be used in the pair correlation function.
Definition: common.h:90
+ Here is the caller graph for this function:

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