Package ffx.potential.nonbonded
Class ParticleMeshEwald
java.lang.Object
ffx.potential.nonbonded.ParticleMeshEwald
- All Implemented Interfaces:
LambdaInterface
This Particle Mesh Ewald class implements PME for the AMOEBA polarizable mutlipole force field in
parallel using a
NeighborList for any Crystal space group. The real space
contribution is contained within this class and the reciprocal space contribution is delegated to
the ReciprocalSpace class.- Author:
- Michael J. Schnieders
derived from:
TINKER code by Jay Ponder, Pengyu Ren and Tom Darden. - See Also:
-
Nested Class Summary
Nested ClassesModifier and TypeClassDescriptionclassThe setFactors(i,k,lambdaMode) method is called every time through the inner PME loops, avoiding an "if (esv)" branch statement.classclass -
Field Summary
FieldsModifier and TypeFieldDescriptionprotected Atom[]An ordered array of atoms in the system.double[][][]Dimensions of [nsymm][xyz][nAtoms].double[][]Direct induced dipoles.double[][]double[][]double[][]doubleCoulomb constant in units of kcal*Ang/(mol*electron^2)double[][][]Fractional multipoles in the global frame with dimensions of [nsymm][nAtoms][10]double[][][]Cartesian multipoles in the global frame with dimensions of [nsymm][nAtoms][10]double[][][]Dimensions of [nsymm][nAtoms][3]double[][][]protected doubleprotected doubleprotected doubleprotected int[][]Polarization groups.protected int[][]protected int[][]The defaults are effectively final, as the implementation of setFactors in the base class is always a no-op.protected intThe number of atoms in the system.int[][][]Neighbor lists, including atoms beyond the real space cutoff.protected doubleprotected doublePermanent multipole energy = permanentRealSpaceEnergy + permanentSelfEnergy + permanentReciprocalEnergy.protected doubleprotected doubleprotected doublePolarization modes include "direct", in which induced dipoles do not interact, and "mutual" that converges the self-consistent field to a tolerance specified by the "polar-eps" keyword.protected doublePolarization energy = inducedRealSpaceEnergy + inducedSelfEnergy + inducedReciprocalEnergy.protected SCFAlgorithmdoubleThe requested permittivity for the solute.protected doubleTotal multipole energy = permanentMultipoleEnergy + polarizationEnergy.double[][]Vacuum induced dipolesdouble[][]double[][][]Vacuum induced dipolesdouble[][][] -
Constructor Summary
ConstructorsConstructorDescriptionParticleMeshEwald(Atom[] atoms, int[] molecule, ForceField forceField, Crystal crystal, NeighborList neighborList, ForceField.ELEC_FORM elecForm, double ewaldCutoff, double gkCutoff, ParallelTeam parallelTeam) ParticleMeshEwald constructor. -
Method Summary
Modifier and TypeMethodDescriptionprotected voidAssignPolarizationGroups.voidattachExtendedSystem(ExtendedSystem system) Attach system with extended variable such as titrations.voiddouble[]computeMoments(Atom[] activeAtoms, boolean forceEnergy) Compute multipole moments for an array of atoms.voiddestroy()doubleenergy(boolean gradient, boolean print) Calculate the PME electrostatic energy.voidReturn the PME AlchemicalParameters.int[][]doubledouble[][][]doubleGet the 2nd partial derivative of the energy with respect to lambda.doublegetdEdL()Get the partial derivative of the energy with respect to lambda.voidgetdEdXdL(double[] gradient) Get the gradient of dEdL with respect to each parameter.doubleReturns the ELEC_FORM.doubledoublegetGK()doublegetGeneralizedKirkwoodEnergy.intgetGKInteractionsdoublegetIndRealEnergy.doublegetIndRecipEnergy.doublegetIndSelfEnergy.intGetter for the fieldinteractions.doubleGet the current value of the state variable.booleanIf true, there are alchemical atoms impacted by the lambda state variable.getMultipoleType(int i) Get the MultipoleType for Atom i.doubledoubledoubledoubledoublegetPermSelfEnergy.doubleint[][]int[][]int[][]doubleGetter for the fieldpolarizationEnergy.getPolarizeType(int i) Get the PolarizeType for Atom i.doubledoublegetGKEnergydoubleGetter for the fieldtotalMultipoleEnergy.voidvoidsetCrystal(Crystal crystal) voidsetLambda(double lambda) Set the current value of the state variable.voidSetter for the fieldpolarization.Methods inherited from class java.lang.Object
clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, waitMethods inherited from interface ffx.potential.bonded.LambdaInterface
dEdLZeroAtEnds
-
Field Details
-
polarization
Polarization modes include "direct", in which induced dipoles do not interact, and "mutual" that converges the self-consistent field to a tolerance specified by the "polar-eps" keyword. -
coordinates
public double[][][] coordinatesDimensions of [nsymm][xyz][nAtoms]. -
neighborLists
public int[][][] neighborListsNeighbor lists, including atoms beyond the real space cutoff. [nsymm][nAtoms][nAllNeighbors] -
globalMultipole
public double[][][] globalMultipoleCartesian multipoles in the global frame with dimensions of [nsymm][nAtoms][10] -
fractionalMultipole
public double[][][] fractionalMultipoleFractional multipoles in the global frame with dimensions of [nsymm][nAtoms][10] -
inducedDipole
public double[][][] inducedDipoleDimensions of [nsymm][nAtoms][3] -
inducedDipoleCR
public double[][][] inducedDipoleCR -
directDipole
public double[][] directDipoleDirect induced dipoles. -
directDipoleCR
public double[][] directDipoleCR -
directField
public double[][] directField -
directFieldCR
public double[][] directFieldCR -
vacuumInducedDipole
public double[][][] vacuumInducedDipoleVacuum induced dipoles -
vacuumInducedDipoleCR
public double[][][] vacuumInducedDipoleCR -
vacuumDirectDipole
public double[][] vacuumDirectDipoleVacuum induced dipoles -
vacuumDirectDipoleCR
public double[][] vacuumDirectDipoleCR -
electric
@FFXProperty(name="electric", propertyGroup=LocalGeometryFunctionalForm, defaultValue="332.063713", description="Specifies a value for the so-called \"electric constant\" allowing conversion unit of electrostatic\npotential energy values from electrons^2/Angstrom to kcal/mol. Internally, FFX stores a default value\nfor this constant as 332.063713 based on CODATA reference values. Since different force fields are\nintended for use with slightly different values, this keyword allows overriding the default value.\n") public double electricCoulomb constant in units of kcal*Ang/(mol*electron^2) -
soluteDielectric
public double soluteDielectricThe requested permittivity for the solute. -
atoms
An ordered array of atoms in the system. -
nAtoms
protected int nAtomsThe number of atoms in the system. -
ip11
protected int[][] ip11Polarization groups. -
ip12
protected int[][] ip12 -
ip13
protected int[][] ip13 -
totalMultipoleEnergy
protected double totalMultipoleEnergyTotal multipole energy = permanentMultipoleEnergy + polarizationEnergy.
This does not include GK. -
permanentMultipoleEnergy
protected double permanentMultipoleEnergyPermanent multipole energy = permanentRealSpaceEnergy + permanentSelfEnergy + permanentReciprocalEnergy. -
permanentRealSpaceEnergy
protected double permanentRealSpaceEnergy -
permanentSelfEnergy
protected double permanentSelfEnergy -
permanentReciprocalEnergy
protected double permanentReciprocalEnergy -
permanentChargeCorrectionEnergy
protected double permanentChargeCorrectionEnergy -
polarizationEnergy
protected double polarizationEnergyPolarization energy = inducedRealSpaceEnergy + inducedSelfEnergy + inducedReciprocalEnergy. -
inducedRealSpaceEnergy
protected double inducedRealSpaceEnergy -
inducedSelfEnergy
protected double inducedSelfEnergy -
inducedReciprocalEnergy
protected double inducedReciprocalEnergy -
scfAlgorithm
-
LambdaDefaults
The defaults are effectively final, as the implementation of setFactors in the base class is always a no-op.
-
-
Constructor Details
-
ParticleMeshEwald
public ParticleMeshEwald(Atom[] atoms, int[] molecule, ForceField forceField, Crystal crystal, NeighborList neighborList, ForceField.ELEC_FORM elecForm, double ewaldCutoff, double gkCutoff, ParallelTeam parallelTeam) ParticleMeshEwald constructor.- Parameters:
atoms- the Atom array to do electrostatics on.molecule- the molecule number for each atom.forceField- the ForceField the defines the electrostatics parameters.crystal- The boundary conditions.neighborList- The NeighborList for both van der Waals and PME.elecForm- The electrostatics functional form.ewaldCutoff- The Ewald real space cutoff.gkCutoff- The generalized Kirkwood cutoff.parallelTeam- A ParallelTeam that delegates parallelization.
-
-
Method Details
-
getElecForm
Returns the ELEC_FORM. -
getLambdaTerm
public boolean getLambdaTerm()If true, there are alchemical atoms impacted by the lambda state variable.- Returns:
- True if there are alchemical atoms impacted by the lambda state variable.
-
getAlchemicalParameters
Return the PME AlchemicalParameters. -
computeInduceDipoleField
public void computeInduceDipoleField() -
destroy
public void destroy() -
energy
public double energy(boolean gradient, boolean print) Calculate the PME electrostatic energy.- Parameters:
gradient- Iftrue, the gradient will be calculated.print- Iftrue, extra logging is enabled.- Returns:
- return the total electrostatic energy (permanent + polarization).
-
expandInducedDipoles
public void expandInducedDipoles() -
getAxisAtoms
public int[][] getAxisAtoms() -
getCavitationEnergy
public double getCavitationEnergy() -
getCoordinates
public double[][][] getCoordinates() -
getDispersionEnergy
public double getDispersionEnergy() -
getEwaldCoefficient
public double getEwaldCoefficient() -
getEwaldCutoff
public double getEwaldCutoff() -
getGK
-
getGKEnergy
public double getGKEnergy()getGeneralizedKirkwoodEnergy.- Returns:
- a double.
-
getGKInteractions
public int getGKInteractions()getGKInteractions- Returns:
- The number of GK interactions.
-
getInteractions
public int getInteractions()Getter for the fieldinteractions.- Returns:
- The number of PME interactions.
-
getLambda
public double getLambda()Get the current value of the state variable.Get the current lambda scale value.
- Specified by:
getLambdain interfaceLambdaInterface- Returns:
- state
-
setLambda
public void setLambda(double lambda) Set the current value of the state variable. May be ignored if lambda is not being applied.Set the electrostatic lambda scaling factor.
- Specified by:
setLambdain interfaceLambdaInterface- Parameters:
lambda- a double.
-
getPermanentEnergy
public double getPermanentEnergy() -
getIndRealEnergy
public double getIndRealEnergy()getIndRealEnergy.- Returns:
- a double.
-
getIndRecipEnergy
public double getIndRecipEnergy()getIndRecipEnergy.- Returns:
- a double.
-
getIndSelfEnergy
public double getIndSelfEnergy()getIndSelfEnergy.- Returns:
- a double.
-
getPermSelfEnergy
public double getPermSelfEnergy()getPermSelfEnergy.- Returns:
- a double.
-
getPermRealEnergy
public double getPermRealEnergy() -
getPermRecipEnergy
public double getPermRecipEnergy() -
getPermanentChargeCorrectionEnergy
public double getPermanentChargeCorrectionEnergy() -
getPolarEps
public double getPolarEps() -
getPolarization11
public int[][] getPolarization11() -
getPolarization12
public int[][] getPolarization12() -
getPolarization13
public int[][] getPolarization13() -
getPolarizationEnergy
public double getPolarizationEnergy()Getter for the fieldpolarizationEnergy.- Returns:
- a double.
-
getPolarizationType
-
getReciprocalSpace
-
getScale14
public double getScale14() -
getSolvationEnergy
public double getSolvationEnergy()getGKEnergy- Returns:
- a double.
-
getMultipoleType
Get the MultipoleType for Atom i.- Parameters:
i- The atom index.- Returns:
- The MultipoleType.
-
getPolarizeType
Get the PolarizeType for Atom i.- Parameters:
i- The atom index.- Returns:
- The PolarizeType.
-
getd2EdL2
public double getd2EdL2()Get the 2nd partial derivative of the energy with respect to lambda.- Specified by:
getd2EdL2in interfaceLambdaInterface- Returns:
- d2EdL2
-
getdEdL
public double getdEdL()Get the partial derivative of the energy with respect to lambda.- Specified by:
getdEdLin interfaceLambdaInterface- Returns:
- dEdL
-
getdEdXdL
public void getdEdXdL(double[] gradient) Get the gradient of dEdL with respect to each parameter.- Specified by:
getdEdXdLin interfaceLambdaInterface- Parameters:
gradient- - A double array of length the number of parameters in the model (commonly 3 * number of atoms).
-
setAtoms
-
setCrystal
-
getTotalMultipoleEnergy
public double getTotalMultipoleEnergy()Getter for the fieldtotalMultipoleEnergy.- Returns:
- a double.
-
setPolarization
Setter for the fieldpolarization.- Parameters:
set- aPolarizationobject.
-
assignPolarizationGroups
protected void assignPolarizationGroups()AssignPolarizationGroups. -
computeMoments
Compute multipole moments for an array of atoms.- Parameters:
activeAtoms- Atom array to consider.forceEnergy- Force calculation of the electrostatic energy (rotate multipoles, perform SCF).
-
attachExtendedSystem
Attach system with extended variable such as titrations.- Parameters:
system- aExtendedSystemobject.
-