iAPBS/NAMD Interface: User Guide
Table of Contents
1 Introduction
The NAMD/APBS interface allows access to APBS functionality from NAMD. The APBS NAMD module can be used to perform implicit solvent MD simulation with NAMD, to write out electrostatic maps for visualization purposes and also to perform MM-PBSA calculations directly with NAMD.
2 Installing the iAPBS Interface
2.1 Requirements
To compile iAPBS/NAMD interface two packages are required:
- APBS (version 1.4 or later).
- NAMD (version 2.11 or later)
2.2 Building iAPBS interface
Building of the iAPBS interface requires compilation and installation of MALOC and APBS libraries. The following describes each step separately (assumes bash as login script, modify appropriately for t/csh).
# create a build directory and cd to it, then: P=`pwd` export APBS_PREFIX=${P} export APBS_SRC=${P}/apbs-pdb2pqr/apbs export MCSH_HOME=/dev/null # get the source git clone https://github.com/Electrostatics/apbs-pdb2pqr cd apbs-pdb2pqr git submodule init git submodule update mkdir apbs/build # configure and build it cd apbs/build cmake -DCMAKE_INSTALL_PREFIX:PATH=${APBS_PREFIX} \ -DCMAKE_BUILD_TYPE=Release -DBUILD_DOC=OFF \ -DBUILD_SHARED_LIBS=OFF -DENABLE_QUIET=ON \ -DENABLE_iAPBS=ON -DBUILD_WRAPPER=ON .. make install
2.3 Building NAMD/APBS module
Currently the NAMD/APBS module is supported on the x86_64 Linux platform only. The Intel compiler is recommended.
cd $APBS_PREFIX #untar NAMD 2.11 source here ln -s NAMD_2.11_Source namd2 export NAMD_SRC=$P/namd2 cd $NAMD_SRC tar xf charm-6.7.0.tar cd charm-6.7.0 ./build charm++ multicore-linux64 --with-production cd $NAMD_SRC wget http://www.ks.uiuc.edu/Research/namd/libraries/fftw-linux-x86_64.tar.gz tar xzf fftw-linux-x86_64.tar.gz mv linux-x86_64 fftw wget http://www.ks.uiuc.edu/Research/namd/libraries/tcl8.5.9-linux-x86_64.tar.gz wget http://www.ks.uiuc.edu/Research/namd/libraries/tcl8.5.9-linux-x86_64-threaded.tar.gz tar xzf tcl8.5.9-linux-x86_64.tar.gz tar xzf tcl8.5.9-linux-x86_64-threaded.tar.gz mv tcl8.5.9-linux-x86_64 tcl mv tcl8.5.9-linux-x86_64-threaded tcl-threaded cd src ln -s $APBS_PREFIX/apbs-pdb2pqr/apbs/contrib/iapbs/modules/NAMD/GlobalMasterAPBS.* . cd ../arch ln -s $APBS_PREFIX/apbs-pdb2pqr/apbs/contrib/iapbs/modules/NAMD/*.apbs . cd $NAMD_SRC patch -p0 < $APBS_PREFIX/apbs-pdb2pqr/apbs/contrib/iapbs/modules/NAMD/patches/namd2-apbs-2.11.patch ./config apbs Linux-x86_64-g++ --charm-arch multicore-linux64 cd Linux-x86_64-g++ make # run tests cd $APBS_PREFIX/apbs-pdb2pqr/apbs/contrib/iapbs/modules/NAMD/test/dipeptide $NAMD_SRC/Linux-x86_64-g++/namd2 apbs-solvation.conf cd ../ubq/ $NAMD_SRC/Linux-x86_64-g++/namd2 ubq.conf
3 Using NAMD/APBS module
The APBS module in NAMD adds implicit solvent simulation environment to the NAMD capabilities. The APBS module can be used for example for implicit solvent minimization and dynamics, calculation and visualization of miscellaneous electrostatic biomolecular properties. Please see the Examples section of this User's Guide.
The following table lists all NAMD/APBS keywords with their description. The left side of the table also lists corresponding APBS keywords which NAMD/APBS keywords mimic very closely. For detailed discussion of APBS keywords please see APBS documentation.
3.1 NAMD/APBS Module Keywords
APBS keyword | NAMD/APBS keyword | Description |
---|---|---|
apbsForces [off] | Turns on the APBS module. | |
apbsPQRFile | PQR file name to be read. | |
mg-auto/mg-para | calc_type [0] | 0: manual MG; 1: autoMG; 2: |
parallel MG | ||
lpbe/nbpe | nonlin [0] | Linear/full Poisson-Boltzmann equation: 0: |
linear; 1: non-linear; 4: size-dependent PBE | ||
bcfl | bcfl [1] | Boundary condition method: 0: zero; 1: sdh; 2: |
mdh; 4: focus | ||
srfm | srfm [1] | Surface calculation method: 0: mol; 1: smol; 2: |
spl2; 3: spl4 | ||
pdie | pdie [2.0] | Solute dielectric |
sdie | sdie [78.4] | Solvent dielectric |
sdens | sdens [10.0] | Vacc sphere density |
srad | srad [1.4] | Solvent radius |
swin | swin [0.3] | Cubic spline window |
temp | temp [298.15] | Temperature (in K) |
gamma | gamma [0.105] | Surface tension for apolar energies/forces |
(in kJ/mol/A2) | ||
chgm | chgm [1] | Charge discretization method: 0: spl0; 1: spl2; 2: |
spl4 | ||
vol | smvolume [10.0] | The parameter smvolume controls the lattice |
size (in Angstroms3) used in the SMPBE formalism. | ||
size | smsize [1000.0] | The parameter smsize controls the relative |
size of the ions (in Angstroms) such that each lattice site can | ||
contain a single ion of volume radius3 or size ions of volume | ||
radius3/size. | ||
calcenergy | calcenergy [1] | Energy calculation flag: 0: Do not |
perform energy calculation; 1: Calculate total energy only; 2: | ||
Calculate per-atom energy components | ||
calcforce | calcforce [0] | Atomic forces calculation: 0: Do not |
perform force calculation; 1: Calculate total force only; 2: | ||
Calculate per-atom force components | ||
calcnpenergy [1] | Calculate nonpolar energy [0, 1]. 0: Don't | |
calculate; 1: SASA-based apolar energy caclulation. | ||
write pot | wpot [off] | Writes electrostatic potential data to |
iapbs-pot.dx in DX format | ||
write charge | wchg [off] | Writes charge data to iapbs-charge.dx in |
DX format | ||
write smol | wsmol [off] | Writes molecular surface data to |
iapbs-smol.dx in DX format | ||
write kappa | wkappa [off] | Writes the ion-accessibility kappa map |
to iapbs-kappa.dx in DX format | ||
write diel | wdiel [off] | Writes dielectric maps to |
iapbs-diel[x,y,z].dx in DX format | ||
read charge | rchg [off] | Reads charge data from iapbs-charge.dx in |
DX format | ||
read kappa | rkappa [off] | Reads the ion-accessibility kappa map |
from iapbs-kappa.dx in DX format | ||
read diel | rdiel [off] | Reads dielectric maps from |
iapbs-diel[x,y,z].dx in DX format | ||
ion | ion charge conc radius | Counterion charges (in e), Counterion |
concentrations (in M), Counterion radii (in A) | ||
dime | dime | Grid dimensions (in x, y and z) |
cmeth | cmeth [1] | Centering method: 0: Center on a point 1: Center |
on a molecule | ||
gcent | center | Grid center if cmeth=0 |
ccmeth | ccmeth [1] | Coarse grid centering method: 0: Center on a |
point 1: Center on a molecule | ||
cgcent | ccenter | Coarse grid center if ccmeth=0 |
fcmeth | fcmeth [1] | Fine grid centering method: 0: Center on a |
point 1: Center on a molecule | ||
fgcent | fcenter | Fine grid center if fcmeth=0 |
grid | grid | Grid spacings |
glen | glen | Grid side lengths |
cglen | cglen | Coarse grid side lengths |
fglen | fglen | Fine grid side lengths |
pdime | pdime | Grid of processors to be used in parallel calculation |
ofrac | ofrac [0.1] | Overlap fraction between processors |
debug [0] | Debugging flag [0-5] | |
verbose [0] | Print (verbosity) flag [0-5] | |
sp_apbs [off] | Perform a single point (a single electrostatic | |
evaluation) APBS calculation | ||
recalculateGrid [off] | Recalculate the grid dimensions on the | |
fly. |
Note
- Values in square brackets [] are defaults.
- When the recalculateGrid and grid keywords are set then the grid size parameters (cglen, fglen and dime) are automatically recalculated on the fly during the simulation. For example if grid is set to [0.5 0.5 0.5] Angstroms the grid size will be automatically adjusted to fit the molecule as the system's dimensions are changed during the simulation (the requested grid spacing 0.5 A will be maintained). This is the recommended option for most of simulations since this prevents the solute to "escape" the pre-set grid.
- To disable creation of io.mc file before attempting any extensive minimization or MD simulation with the APBS module please set the
MCSH_HOME
environment variable to/dev/null
(export MCSH_HOME=/dev/null
).
4 Examples of using the NAMD/APBS module
The APBS module is initialized with keyword apbsForces
. The APBS
calculation set up is specified in the apbsForcesConfig
section. Please note that charges and radii definition for the
electrostatic calculation must be read from a PQR file. See APBS
documentation for PQR file format and description. APBS distribution
also contains tools for generating valid PQR files (from a PDB file,
for example). The order of atoms in the PQR file must be the same as
in the associated .top file. For list of APBS-related keywords see
NAMD/APBS keywords table.
4.1 Solvation energy calculation
This is an example of single point solvation energy calculation. The electrostatic calculation is done on 0.53 A numerical grid. The external charges and radii are read from dipeptide.pqr file. The final solvation energy is printed in the "MISC" energy column (as a sum of electrostatic and non-polar energies).
amber on parmfile dipeptide.top ambercoor dipeptide.crd temperature 300 exclude scaled1-4 1-4scaling 0.8333333 switching on switchDist 9 cutoff 10 pairListDist 11 outputname output outputEnergies 1 outputTiming 100 dcdFreq 500 restartFreq 500 wrapWater on wrapNearest on langevin on langevinDamping 2 langevinHydrogen no langevinTemp 300 apbsForces on apbsPQRFile dipeptide.pqr apbsForcesConfig { calc_type 0 # mg-manual grid 0.5 0.5 0.5 recalculateGrid on srfm 2 chgm 1 bcfl 1 debug 1 verbose 5 pdie 2.0000 sdie 78.5400 sdens 10.00 srad 1.40 swin 0.30 temp 298.15 gamma 0.105 sp_apbs off wpot off } minimize 0
4.2 Molecular dynamics in implicit solvent using APBS
This calculation performs a MD simulation in implicit solvent (water in this case). The charges and radii definition is read from an external file (dipeptide.pqr) and the electrostatic calculation is done on 333 points of numerical grid with the grid dimensions specified in the input file.
amber on parmfile dipeptide.top ambercoor dipeptide.crd temperature 300 exclude scaled1-4 1-4scaling 0.8333333 switching on switchDist 9 cutoff 10 pairListDist 11 outputname output outputEnergies 1 outputTiming 100 dcdFreq 500 restartFreq 500 wrapWater on wrapNearest on langevin off langevinDamping 2 langevinHydrogen no langevinTemp 300 apbsForces on apbsPQRFile dipeptide.pqr apbsForcesConfig { dime 33 33 33 cglen 17.0071 13.8706 12.3012 fglen 17.0071 13.8706 12.3012 srfm 2 chgm 1 bcfl 1 debug 0 pdie 2.0000 sdie 78.5400 sdens 10.00 srad 1.40 swin 0.30 temp 298.15 gamma 0.105 sp_apbs off wpot off } numsteps 100
4.3 Calculation and visualization of electrostatic potential
This calculation writes out the calculated electrostatic potential to
a file (iapbs-pot.dx
). This potential can be visualized using vmd and
pymol. The grid dimensions are 1.03 A.
amber on parmfile dipeptide.top ambercoor dipeptide.crd temperature 300 exclude scaled1-4 1-4scaling 0.8333333 switching on switchDist 9 cutoff 10 pairListDist 11 outputname output outputEnergies 1 outputTiming 100 dcdFreq 500 restartFreq 500 wrapWater on wrapNearest on langevin on langevinDamping 2 langevinHydrogen no langevinTemp 300 apbsForces on apbsPQRFile dipeptide.pqr apbsForcesConfig { dime 0 0 0 grid 1.0 1.0 1.0 srfm 2 chgm 1 bcfl 1 debug 0 verbose 2 pdie 2.0000 sdie 78.5400 sdens 10.00 srad 1.40 swin 0.30 temp 298.15 gamma 0.105 sp_apbs on wpot on } minimize 0
Note
To determine the correct size of the numerical grid enveloping the studied molecule (parameters
cglen
andfglen
) you can usepsize.py
tool from the APBS distribution.