EEFSOLV
This is part of the colvar module

Calculates EEF1 solvation free energy for a group of atoms.

EEF1 is a solvent-accessible surface area based model, where the free energy of solvation is computed using a pairwise interaction term for non-hydrogen atoms:

\[ \Delta G^\mathrm{solv}_i = \Delta G^\mathrm{ref}_i - \sum_{j \neq i} f_i(r_{ij}) V_j \]

where \(\Delta G^\mathrm{solv}_i\) is the free energy of solvation, \(\Delta G^\mathrm{ref}_i\) is the reference solvation free energy, \(V_j\) is the volume of atom \(j\) and

\[ f_i(r) 4\pi r^2 = \frac{2}{\sqrt{\pi}} \frac{\Delta G^\mathrm{free}_i}{\lambda_i} \exp\left\{ - \frac{(r-R_i)^2}{\lambda^2_i}\right\} \]

where \(\Delta G^\mathrm{free}_i\) is the solvation free energy of the isolated group, \(\lambda_i\) is the correlation length equal to the width of the first solvation shell and \(R_i\) is the van der Waals radius of atom \(i\).

The output from this collective variable, the free energy of solvation, can be used with the BIASVALUE keyword to provide implicit solvation to a system. All parameters are designed to be used with a modified CHARMM36 force field. It takes only non-hydrogen atoms as input, these can be conveniently specified using the GROUP action with the NDX_GROUP parameter. To speed up the calculation, EEFSOLV internally uses a neighbor list with a cutoff dependent on the type of atom (maximum of 1.95 nm). This cutoff can be extended further by using the NL_BUFFER keyword.

Examples
Click on the labels of the actions for more information on what each action computes
tested on v2.9
#SETTINGS MOLFILE=regtest/basic/rt77/peptide.pdb
MOLINFO MOLTYPEcompulsory keyword ( default=protein )
what kind of molecule is contained in the pdb file - usually not needed since protein/RNA/DNA
are compatible =protein STRUCTUREcompulsory keyword 
a file in pdb format containing a reference structure. =peptide.pdb 
WHOLEMOLECULES ENTITY0the atoms that make up a molecule that you wish to align. =1-111 
# This allows us to select only non-hydrogen atoms
#SETTINGS AUXFILE=regtest/basic/rt77/index.ndx
protein-h: GROUP NDX_FILEthe name of index file (gromacs syntax) =index.ndx NDX_GROUPthe name of the group to be imported (gromacs syntax) - first group found is used
by default =Protein-H 
# We extend the cutoff by 0.1 nm and update the neighbor list every 40 steps
solv: EEFSOLV ATOMSThe atoms to be included in the calculation, e.g. =protein-h 
# Here we actually add our calculated energy back to the potential
bias: BIASVALUE ARGthe input for this action is the scalar output from one or more other actions. =solv 
PRINT ARGthe input for this action is the scalar output from one or more other actions. =solv FILEthe name of the file on which to output these quantities =SOLV 
Glossary of keywords and components
The atoms involved can be specified using
ATOMS The atoms to be included in the calculation, e.g. the whole protein.. For more information on how to specify lists of atoms see Groups and Virtual Atoms
Compulsory keywords
NL_BUFFER ( default=0.1 ) The buffer to the intrinsic cutoff used when calculating pairwise interactions.
NL_STRIDE ( default=40 ) The frequency with which the neighbor list is updated.
Options
NUMERICAL_DERIVATIVES ( default=off ) calculate the derivatives for these quantities numerically
NOPBC ( default=off ) ignore the periodic boundary conditions when calculating distances
SERIAL ( default=off ) Perform the calculation in serial - for debug purpose
TEMP_CORRECTION

( default=off ) Correct free energy of solvation constants for temperatures different from 298.15 K