espressopp.integrator.PIAdressIntegrator¶
The PIAdressIntegrator implements the integration method for Hamiltonian Adaptive Resolution Path Integral Simulations proposed in J. Chem. Phys 147, 244104 (2017) (PI-AdResS). It can be used to run path integral molecular dynamics as well as ring polymer and centroid molecular dynamics in a quantum-classical adaptive resolution fashion, using different empirical force fields. To facilitate an efficient integration, the integrator uses a 3-layer RESPA (J. Chem. Phys. 97, 1990 (1992)) multiple timestepping scheme (inner level: intraatomic spring forces between the Trotter beads. medium level: interatomic bonded forces. outer level: interatomic non-bonded forces). Importantly, the integrator should only be used in combination with PI-AdResS interactions. Furthermore, the integrator has its own thermostat (Langevin), and the only extensions that should be used with it are the Free Energy Compensation (FreeEnergyCompensation) and the Thermodynamic Force (TDforce).
Example:
>>> integrator = espressopp.integrator.PIAdressIntegrator(system, verletlist, timestep_short, timesteps_outerlevel, timesteps_centrallevel, nTrotter, realkinmass, constkinmass, temperature, gamma, centroidthermostat, CMDparameter, PILE, PILElambda, clmassmultiplier, speedupFreezeRings, KTI)
>>> ...
>>> integrator.run(nsteps)
-
espressopp.integrator.
PIAdressIntegrator
(system, verletlist, timestep, sSteps, mSteps, nTrotter, realKinMass, constKinMass, temperature, gamma, centroidThermostat, CMDparameter, PILE, PILElambda, CLmassmultiplier, speedup, KTI)¶ Constructs the PIAdressIntegrator object. Note that all parameters can also be set and fetched via setter and getter functions. Additionally, all parameters except the system and the Verletlist are implemented as class variables that can be directly accessed and modified.
Parameters: - system (shared_ptr<System>) – system object
- verletlist (shared_ptr<VerletListAdress>) – Verletlist object. Should be an AdResS Verletlist
- timestep (real) – (default: 0.0) the inner (shortest) timestep for the calculation of the intraatomic spring forces between the Trotter beads
- sSteps (int) – (default: 1) multiplier to construct medium timestep (interatomic bonded forces) as mediumstep = sSteps * timestep
- mSteps (int) – (default: 1) multiplier to construct longest timestep (interatomic non-bonded forces) as longstep = mSteps * sSteps * timestep
- nTrotter (int) – (default: 32) Trotter number. Should be even and greather than zero.
- realKinMass (bool) – (default: True) Flag to choose whether to use real kinetic masses. If False, the higher modes’ kinetic masses are multiplied with their corresponding eigenvalues of the normal mode transformation. In this way, all higher modes oscillate with the same frequency. If True, we use the kinetic masses for the higher modes which corresponding to the real dynamics (see J. Chem. Phys 147, 244104 (2017) for details)
- constKinMass (bool) – (default: False) If False, the higher modes’ kinetic masses also adaptively change (AKM scheme in J. Chem. Phys 147, 244104 (2017)). If True, the higher modes’ kinetic masses are constant throughout the system (CKM scheme in J. Chem. Phys 147, 244104 (2017))
- temperature (real) – (default: 2.494353 - this corresponds to 300 Kelvin) the temperature in gromacs units (Boltzmann constant kb is 1)
- gamma (real) – (default: 1.0) the Langevin thermostat’s friction parameter in 1/ps
- centroidThermostat (bool) – (default: True) If True, the centroid mode is also thermostated, otherwise only the higher modes’ (relevant for centroid molecular dynamics)
- CMDparameter (real) – (default: 1.0) The gamma^2 parameter used in centroid molecular dynamics. The higher modes’ kinetic masses are rescaled by CMDparameter
- PILE (bool) – (default: True) If True, the higher modes are thermostated according to the PILE scheme by Ceriotti et al. (J. Chem. Phys 133, 124104 (2010)). Only makes sense in combination when using real kinetic masses (realKinMass = True)
- PILElambda (real) – (default: 0.5) lambda parameter to rescale the friction matrix. Default should be good for most applications (J. Chem. Phys 140, 234116 (2014))
- CLmassmultiplier (real) – (default: 100.0) multiplier by which the higher modes’ spring masses (if constKinMass = False also the kinetic masses) are increased in the classical region
- speedup (bool) – (default: True) If True, the higher modes’ are not integrated in the classical region and also the intraatomistic forces between the Trotter beads are not calculated in the classical region
- KTI (bool) – (default: False) If True, the particles’ resolution parameters and adaptive masses are not updated but can be set by hand everywhere. This is necessary when running Kirkwood Thermodynamic Integration (KTI)
-
espressopp.integrator.PIAdressIntegrator.
setVerletList
(verletlist)¶ Sets the VerletList.
Parameters: verletlist (espressopp.VerletListAdress) – The VerletListAdress object.
-
espressopp.integrator.PIAdressIntegrator.
getVerletList
()¶ Gets the VerletList.
Returns: the Adress VerletList Return type: shared_ptr<VerletListAdress>
-
espressopp.integrator.PIAdressIntegrator.
setTimeStep
(timestep)¶ Sets the inner (shortest) timestep.
Parameters: timestep (real) – the inner timestep
-
espressopp.integrator.PIAdressIntegrator.
getTimeStep
()¶ Gets the inner (shortest) timestep.
Returns: the inner timestep Return type: real
-
espressopp.integrator.PIAdressIntegrator.
setsStep
(sSteps)¶ Sets the multiplier to construct medium timestep (interatomic bonded forces) as mediumstep = sSteps * timestep.
Parameters: sSteps (int) – multiplier to construct medium timestep
-
espressopp.integrator.PIAdressIntegrator.
getsStep
()¶ Gets the multiplier to construct medium timestep (interatomic bonded forces) as mediumstep = sSteps * timestep.
Returns: multiplier to construct medium timestep Return type: int
-
espressopp.integrator.PIAdressIntegrator.
setmStep
(mSteps)¶ Sets the multiplier to construct longest timestep (interatomic non-bonded forces) as longstep = mSteps * sSteps * timestep.
Parameters: mSteps (int) – multiplier to construct longest timestep
-
espressopp.integrator.PIAdressIntegrator.
getmStep
()¶ Gets the multiplier to construct longest timestep (interatomic non-bonded forces) as longstep = mSteps * sSteps * timestep.
Returns: multiplier to construct longest timestep Return type: int
-
espressopp.integrator.PIAdressIntegrator.
setNtrotter
(nTrotter)¶ Sets the Trotter number nTrotter. Should be even and greather than zero. Note that when calling this function, also the normal mode transformation matrix and the eigenvalues are recalculated.
Parameters: ntrotter (int) – the Trotter number
-
espressopp.integrator.PIAdressIntegrator.
getNtrotter
()¶ Gets the Trotter number nTrotter.
Returns: the Trotter number Return type: int
-
espressopp.integrator.PIAdressIntegrator.
setRealKinMass
(realKinMass)¶ Sets the real kinetic mass flag.
Parameters: realKinMass (bool) – the real kinetic mass flag
-
espressopp.integrator.PIAdressIntegrator.
getRealKinMass
()¶ Gets the real kinetic mass flag.
Returns: the real kinetic mass flag Return type: bool
-
espressopp.integrator.PIAdressIntegrator.
setConstKinMass
(constKinMass)¶ Sets the constant kinetic mass flag.
Parameters: constKinMass (bool) – the constant kinetic mass flag
-
espressopp.integrator.PIAdressIntegrator.
getConstKinMass
()¶ Gets the constant kinetic mass flag.
Returns: the constant kinetic mass flag Return type: bool
-
espressopp.integrator.PIAdressIntegrator.
setTemperature
(temperature)¶ Sets the temperature (gromacs units with kb = 1).
Parameters: temperature (real) – the temperature
-
espressopp.integrator.PIAdressIntegrator.
getTemperature
()¶ Gets the temperature (gromacs units with kb = 1).
Returns: the temperature Return type: real
-
espressopp.integrator.PIAdressIntegrator.
setGamma
(gamma)¶ Sets the friction constant gamma (in 1/ps).
Parameters: gamma (real) – the friction constant gamma
-
espressopp.integrator.PIAdressIntegrator.
getGamma
()¶ Gets the friction constant gamma (in 1/ps).
Returns: the friction constant gamma Return type: real
-
espressopp.integrator.PIAdressIntegrator.
setCentroidThermostat
(centroidThermostat)¶ Sets the centroid thermostat flag.
Parameters: centroidThermostat (bool) – the centroid thermostat flag
-
espressopp.integrator.PIAdressIntegrator.
getCentroidThermostat
()¶ Gets the centroid thermostat flag.
Returns: the centroid thermostat flag Return type: bool
-
espressopp.integrator.PIAdressIntegrator.
setCMDparameter
(CMDparameter)¶ Sets the centroid molecular dynamics parameter gamma^2 for scaling the kinetic mass.
Parameters: CMDparameter (real) – the CMD parameter gamma^2
-
espressopp.integrator.PIAdressIntegrator.
getCMDparameter
()¶ Gets the centroid molecular dynamics parameter gamma^2 for scaling the kinetic mass.
Returns: the CMD parameter gamma^2 Return type: real
-
espressopp.integrator.PIAdressIntegrator.
setPILE
(PILE)¶ Sets the PILE flag.
Parameters: PILE (bool) – the PILE flag
-
espressopp.integrator.PIAdressIntegrator.
getPILE
()¶ Gets the PILE flag.
Returns: the PILE flag Return type: bool
-
espressopp.integrator.PIAdressIntegrator.
setPILElambda
(PILElambda)¶ Sets the scaling parameter lambda of the PILE thermostat.
Parameters: PILElambda (real) – the scaling parameter lambda
-
espressopp.integrator.PIAdressIntegrator.
getPILElambda
()¶ Gets the scaling parameter lambda of the PILE thermostat.
Returns: the scaling parameter lambda Return type: real
-
espressopp.integrator.PIAdressIntegrator.
setClmassmultiplier
(CLmassmultiplier)¶ Sets the multiplier for the higher modes’ spring masses in the classical region.
Parameters: CLmassmultiplier (real) – the classical spring mass multiplier
-
espressopp.integrator.PIAdressIntegrator.
getClmassmultiplier
()¶ Gets the multiplier for the higher modes’ spring masses in the classical region.
Returns: the classical spring mass multiplier Return type: real
-
espressopp.integrator.PIAdressIntegrator.
setSpeedup
(speedup)¶ Sets the speedup flag.
Parameters: speedup (bool) – the speedup flag
-
espressopp.integrator.PIAdressIntegrator.
getSpeedup
()¶ Gets the speedup flag.
Returns: the speedup flag Return type: bool
-
espressopp.integrator.PIAdressIntegrator.
setKTI
(KTI)¶ Sets the KTI flag.
Parameters: speedup (bool) – the KTI flag
-
espressopp.integrator.PIAdressIntegrator.
getKTI
()¶ Gets the KTI flag.
Returns: the KTI flag Return type: bool
-
espressopp.integrator.PIAdressIntegrator.
getVerletlistBuilds
()¶ Gets the number of Verletlist builds.
Returns: number of Verletlist builds Return type: int
-
espressopp.integrator.PIAdressIntegrator.
computeRingEnergy
()¶ Calculates the total configurational energy of all ring polymers in the system based on the springs between the Trotter beads (calculation done using mode coordinates).
Returns: total configurational ring polymer energy Return type: real
-
espressopp.integrator.PIAdressIntegrator.
computeRingEnergyRaw
()¶ Calculates the total configurational energy of all ring polymers in the system based on the springs between the Trotter beads (calculation done using the Trotter beads’ real space positions).
Returns: total configurational ring polymer energy Return type: real
-
espressopp.integrator.PIAdressIntegrator.
computeKineticEnergy
()¶ Calculates the total kinetic energy using the modes’ momenta.
Returns: total kinetic energy Return type: real
-
espressopp.integrator.PIAdressIntegrator.
computePositionDrift
(parttype)¶ Calculates the average drift force due to the position-dependent spring masses (see Section 5.C. Eq. 63 in J. Chem. Phys 147, 244104 (2017)) on particles of type parttype. To be used during KTI for construction of free energy compensation.
Parameters: parttype (int) – the particle or atom type Returns: average drift force due to the position-dependent spring masses Return type: real
-
espressopp.integrator.PIAdressIntegrator.
computeMomentumDrift
(parttype)¶ Calculates the average drift force due to the position-dependent kinetic masses (see Section 5.C. Eq. 62 in J. Chem. Phys 147, 244104 (2017)) on particles of type parttype. To be used during KTI for construction of free energy compensation.
Parameters: parttype (int) – the particle or atom type Returns: average drift force due to the position-dependent kinetic masses Return type: real
-
class
espressopp.integrator.PIAdressIntegrator.
PIAdressIntegratorLocal
(system, verletlist, timestep=0.0, sSteps=1, mSteps=1, nTrotter=32, realKinMass=True, constKinMass=False, temperature=2.494353, gamma=1.0, centroidThermostat=True, CMDparameter=1.0, PILE=True, PILElambda=0.5, CLmassmultiplier=100.0, speedup=True, KTI=False)¶ The (local) PIAdress Integrator.