This is part of the sasa module | |
It is only available if you configure PLUMED with ./configure –enable-modules=sasa . Furthermore, this feature is still being developed so take care when using it and report any problems on the mailing list. |
Calculates the solvent accessible surface area (SASA) of a protein molecule, or other properties related to it.
The atoms for which the SASA is desired should be indicated with the keyword ATOMS, and a pdb file of the protein must be provided in input with the MOLINFO keyword. The algorithm described in [65] is used for the calculation. The radius of the solvent is assumed to be 0.14 nm, which is the radius of water molecules. Using the keyword NL_STRIDE it is also possible to specify the frequency with which the neighbor list for the calculation of SASA is updated (the default is every 10 steps).
Different properties can be calculated and selected using the TYPE keyword:
1) the total SASA (TOTAL);
2) the free energy of transfer for the protein according to the transfer model (TRANSFER). This keyword can be used, for instance, to compute the transfer of a protein to different temperatures, as detailed in [6], or to different pressures, as detailed in [5], or to different osmolyte solutions, as detailed in [7].
When the TRANSFER keyword is used, a file with the free energy of transfer values for the sidechains and backbone atoms should be provided (using the keyword DELTAGFILE). Such file should have the following format:
----------------Sample DeltaG.dat file--------------------- ALA 0.711019999999962 ARG -2.24832799999996 ASN -2.74838799999999 ASP -2.5626376 CYS 3.89864000000006 GLN -1.76192 GLU -2.38664400000002 GLY 0 HIS -3.58152799999999 ILE 2.42634399999986 LEU 1.77233599999988 LYS -1.92576400000002 MET -0.262827999999956 PHE 1.62028800000007 PRO -2.15598800000001 SER -1.60934800000004 THR -0.591559999999987 TRP 1.22936000000027 TYR 0.775547999999958 VAL 2.12779200000011 BACKBONE 1.00066920000002 -----------------------------------------------------------
where the second column is the free energy of transfer for each sidechain/backbone, in kJ/mol.
A Python script for the computation of free energy of transfer values to describe the effect of osmolyte concentration, temperature and pressure (according to [7], [6] and [5]) is freely available at https://github.com/andrea-arsiccio/DeltaG-calculation. The script automatically outputs a DeltaG.dat file compatible with this SASA module.
If the DELTAGFILE is not provided, the program computes the free energy of transfer values as if they had to take into account the effect of temperature according to approaches 2 or 3 in the paper [6]. Please read and cite this paper if using the transfer model for computing the effect of temperature in implicit solvent simulations. For this purpose, the keyword APPROACH should be added, and set to either 2 or 3.
The SASA usually makes sense when atoms used for the calculation are all part of the same molecule. When running with periodic boundary conditions, the atoms should be in the proper periodic image. This is done automatically since PLUMED 2.2, by considering the ordered list of atoms and rebuilding the broken entities using a procedure that is equivalent to that done in WHOLEMOLECULES. Notice that rebuilding is local to this action. This is different from WHOLEMOLECULES which actually modifies the coordinates stored in PLUMED.
In case you want to recover the old behavior you should use the NOPBC flag. In that case you need to take care that atoms are in the correct periodic image.
The SASA may also be computed using the SASA_LCPO collective variable, which makes use of the LCPO algorithm [143]. SASA_LCPO is more accurate then SASA_HASEL, but the computation is slower.
The following input tells plumed to print the total SASA for atoms 10 to 20 in a protein chain.
sasa: SASA_HASELTYPE=TOTALcompulsory keyword ( default=TOTAL ) The type of calculation you want to perform.ATOMS=10-20the group of atoms that you are calculating the SASA for.NL_STRIDE=10 PRINTcompulsory keyword The frequency with which the neighbor list for the calculation of SASA is updated.ARG=sasathe input for this action is the scalar output from one or more other actions.STRIDE=1compulsory keyword ( default=1 ) the frequency with which the quantities of interest should be outputFILE=colvarthe name of the file on which to output these quantities
The following input tells plumed to compute the transfer free energy for the protein chain containing atoms 10 to 20. Such transfer free energy is then used as a bias in the simulation (e.g., implicit solvent simulations). The free energy of transfer values are read from a file called DeltaG.dat.
sasa: SASA_HASELTYPE=TRANSFERcompulsory keyword ( default=TOTAL ) The type of calculation you want to perform.ATOMS=10-20the group of atoms that you are calculating the SASA for.NL_STRIDE=10compulsory keyword The frequency with which the neighbor list for the calculation of SASA is updated.DELTAGFILE=DeltaG.dat bias: BIASVALUEa file containing the free energy of transfer values for backbone and sidechains atoms.ARG=sasa PRINTthe input for this action is the scalar output from one or more other actions.ARG=sasa,bias.*the input for this action is the scalar output from one or more other actions.STRIDE=1compulsory keyword ( default=1 ) the frequency with which the quantities of interest should be outputFILE=colvarthe name of the file on which to output these quantities
The following input tells plumed to compute the transfer free energy for the protein chain containing atoms 10 to 20. Such transfer free energy is then used as a bias in the simulation (e.g., implicit solvent simulations). The free energy of transfer values are computed according to [6], and take into account the effect of temperature using approach 2 as described in the paper.
sasa: SASA_HASELTYPE=TRANSFERcompulsory keyword ( default=TOTAL ) The type of calculation you want to perform.ATOMS=10-20the group of atoms that you are calculating the SASA for.NL_STRIDE=10compulsory keyword The frequency with which the neighbor list for the calculation of SASA is updated.APPROACH=2 bias: BIASVALUEeither approach 2 or 3.ARG=sasa PRINTthe input for this action is the scalar output from one or more other actions.ARG=sasa,bias.*the input for this action is the scalar output from one or more other actions.STRIDE=1compulsory keyword ( default=1 ) the frequency with which the quantities of interest should be outputFILE=colvarthe name of the file on which to output these quantities
ATOMS | the group of atoms that you are calculating the SASA for. For more information on how to specify lists of atoms see Groups and Virtual Atoms |
TYPE | ( default=TOTAL ) The type of calculation you want to perform. Can be TOTAL or TRANSFER |
NL_STRIDE | The frequency with which the neighbor list for the calculation of SASA is updated. |
NUMERICAL_DERIVATIVES | ( default=off ) calculate the derivatives for these quantities numerically |
NOPBC | ( default=off ) ignore the periodic boundary conditions when calculating distances |
DELTAGFILE | a file containing the free energy of transfer values for backbone and sidechains atoms. Necessary only if TYPE = TRANSFER. A Python script for the computation of free energy of transfer values to describe the effect of osmolyte concentration, temperature and pressure is freely available at https://github.com/andrea-arsiccio/DeltaG-calculation. The script automatically outputs a DeltaG.dat file compatible with this SASA module. If TYPE = TRANSFER and no DELTAGFILE is provided, the free energy values are those describing the effect of temperature, and are computed using the temperature value passed by the MD engine |
APPROACH | either approach 2 or 3. Necessary only if TYPE = TRANSFER and no DELTAGFILE is provided. If TYPE = TRANSFER and no DELTAGFILE is provided, the free energy values are those describing the effect of temperature, and the program must know if approach 2 or 3 (as described in Arsiccio and Shea, Protein Cold Denaturation in Implicit Solvent Simulations: A Transfer Free Energy Approach, J. Phys. Chem. B, 2021) needs to be used to compute them |