Module |
secondarystructure |
Description |
Usage |
Probe the alpha helical content of a protein structure. |
  |
output value |
type |
if LESS_THAN is present the RMSD distance between each residue and the ideal alpha helix |
scalar/vector |
Output components
This action calculates the values in the following table. These values can be referenced elsewhere in the input by using this Action's label followed by a dot and the name of the value required from the list below.
Name |
Type |
Description |
struct |
vector |
the vectors containing the rmsd distances between the residues and each of the reference structures |
lessthan |
scalar |
the number blocks of residues that have an RMSD from the secondary structure that is less than the threshold |
The atoms that serve as the input for this action are specified using one or more of the keywords in the following table.
Keyword |
Type |
Description |
RESIDUES |
atoms |
this command is used to specify the set of residues that could conceivably form part of the secondary structure |
Further details and examples
Probe the alpha helical content of a protein structure.
Any chain of six contiguous residues in a protein chain can form an alpha helix. This
colvar thus generates the set of all possible six residue sections and calculates
the DRMSD or RMSD distance between the configuration in which the residues find themselves
and an idealized alpha helical structure. These distances can be calculated by either
aligning the instantaneous structure with the reference structure and measuring each
atomic displacement or by calculating differences between the set of inter-atomic
distances in the reference and instantaneous structures.
In the original paper on this method (see below) the authors used the set of distances from the alpha helix configurations to measure
the number of segments that have an alpha helical configuration. This is done by calculating
the following sum of functions of the rmsd distances:
s=∑i1−(ri−d0r0)n1−(ri−d0r0)m
where the sum runs over all possible segments of alpha helix. By default the
NN, MM and D_0 parameters are set equal to those used in the paper cited below. The R_0
parameter must be set by the user - the value used in the paper cited below was 0.08 nm.
If you change the function in the above sum you can calculate quantities such as the average
distance from a purely the alpha helical configuration or the distance between the set of
residues that is closest to an alpha helix and the reference configuration. To do these sorts of
calculations you can use the AVERAGE and MIN keywords. In addition you can use the LESS_THAN
keyword if you would like to change the form of the switching function. If you use any of these
options you no longer need to specify NN, R_0, MM and D_0.
The following input calculates the number of six residue segments of
protein that are in an alpha helical configuration.
Click on the labels of the actions for more information on what each action computes
MOLINFOThis command is used to provide information on the molecules that are present in your system. More details STRUCTUREa file in pdb format containing a reference structure= ATOM 1 HH31 ACE 1 -9.105 -2.402 21.804 0.00 0.00
ATOM 2 CH3 ACE 1 -8.930 -3.352 22.308 0.00 0.00
ATOM 3 HH32 ACE 1 -9.504 -3.501 23.223 0.00 0.00
ATOM 4 HH33 ACE 1 -9.067 -4.173 21.604 0.00 0.00
ATOM 5 C ACE 1 -7.450 -3.303 22.659 0.00 0.00
...
ATOM 129 CH3 NME 14 -1.690 -7.089 20.487 0.00 0.00
ATOM 130 HH31 NME 14 -1.873 -6.629 21.459 0.00 0.00
ATOM 131 HH32 NME 14 -1.152 -6.284 19.986 0.00 0.00
ATOM 132 HH33 NME 14 -1.135 -8.020 20.602 0.00 0.00
END
The MOLINFO action with label calculates somethingalphaThe ALPHARMSD action with label alpha calculates the following quantities: Quantity | Type | Description |
alpha | scalar | if LESS_THAN is present the RMSD distance between each residue and the ideal alpha helix. If LESS_THAN is not present the number of residue segments where the structure is similar to an alpha helix |
: ALPHARMSDProbe the alpha helical content of a protein structure. This action is a shortcut and it has hidden defaults. More details RESIDUESthis command is used to specify the set of residues that could conceivably form part of the secondary structure=all R_0The r_0 parameter of the switching function=0.1
alpha: ALPHARMSDProbe the alpha helical content of a protein structure. This action is a shortcut and uses the defaults shown here. More details RESIDUESthis command is used to specify the set of residues that could conceivably form part of the secondary structure=all R_0The r_0 parameter of the switching function=0.1 NN The n parameter of the switching function=8 D_0 The d_0 parameter of the switching function=0.0 MM The m parameter of the switching function=12 TYPE the manner in which RMSD alignment is performed=DRMSD
# alpha: ALPHARMSD RESIDUES=all R_0=0.1
alpha_rmsdThe SECONDARY_STRUCTURE_RMSD action with label alpha_rmsd calculates the following quantities: Quantity | Type | Description |
alpha_rmsd | vector | a vector containing the rmsd distance between each of the residue segments and the reference structure |
: SECONDARY_STRUCTURE_RMSDCalclulate the distance between segments of a protein and a reference structure of interest More details BONDLENGTHthe length to use for bonds=0.17 ATOMSthis is the full list of atoms that we are investigating=7,9,11,15,16,17,19,21,25,26,27,29,31,35,36,37,39,41,45,46,47,49,51,55,56,57,59,61,65,66,67,69,71,75,76,77,79,81,85,86,87,89,91,95,96,97,99,101,105,106,107,109,111,115,116,117,119,121,125,126 TYPE the manner in which RMSD alignment is performed=DRMSD SEGMENT1this is the lists of atoms in the segment that are being considered=0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,17,18,19,20,21,22,23,24,25,26,27,28,29 SEGMENT2this is the lists of atoms in the segment that are being considered=5,6,7,8,9,10,11,12,13,14,15,16,17,18,19,20,21,22,23,24,25,26,27,28,29,30,31,32,33,34 SEGMENT3this is the lists of atoms in the segment that are being considered=10,11,12,13,14,15,16,17,18,19,20,21,22,23,24,25,26,27,28,29,30,31,32,33,34,35,36,37,38,39 SEGMENT4this is the lists of atoms in the segment that are being considered=15,16,17,18,19,20,21,22,23,24,25,26,27,28,29,30,31,32,33,34,35,36,37,38,39,40,41,42,43,44 SEGMENT5this is the lists of atoms in the segment that are being considered=20,21,22,23,24,25,26,27,28,29,30,31,32,33,34,35,36,37,38,39,40,41,42,43,44,45,46,47,48,49 STRUCTURE1the reference structure=0.733,0.519,5.298,1.763,0.81,4.301,3.166,0.543,4.881,1.527,-0.045,3.053,1.646,0.436,1.928,1.18,-1.312,3.254,0.924,-2.203,2.126,0.65,-3.626,2.626,-0.239,-1.711,1.261,-0.19,-1.815,0.032,-1.28,-1.172,1.891,-2.416,-0.661,1.127,-3.548,-0.217,2.056,-1.964,0.529,0.276,-2.364,0.659,-0.88,-1.13,1.391,0.856,-0.62,2.565,0.148,0.228,3.439,1.077,0.231,2.129,-1.032,0.179,2.733,-2.099,1.028,1.084,-0.833,1.872,0.593,-1.919,2.85,-0.462,-1.397,1.02,0.02,-3.049,1.317,0.227,-4.224,-0.051,-0.684,-2.696,-0.927,-1.261,-3.713,-1.933,-2.219,-3.074,-1.663,-0.171,-4.475,-1.916,-0.296,-5.673
alpha_ltThe LESS_THAN action with label alpha_lt calculates the following quantities: Quantity | Type | Description |
alpha_lt | vector | the vector obtained by doing an element-wise application of a function that is one if the input is less than a threshold to the input vectors |
: LESS_THANUse a switching function to determine how many of the input variables are less than a certain cutoff. More details ARGthe values input to this function=alpha_rmsd SWITCHThis keyword is used if you want to employ an alternative to the continuous swiching function defined above={RATIONAL R_0=0.1 D_0=0.0 NN=8 MM=12}
alphaThe SUM action with label alpha calculates the following quantities: Quantity | Type | Description |
alpha | scalar | the sum of all the elements in the input vector |
: SUMCalculate the sum of the arguments More details ARGthe values input to this function=alpha_lt PERIODICif the output of your function is periodic then you should specify the periodicity of the function=NO
# --- End of included input --- PRINTPrint quantities to a file. More details ARGthe labels of the values that you would like to print to the file=alpha FILEthe name of the file on which to output these quantities=colvar
Here the same is done use RMSD instead of DRMSD
Click on the labels of the actions for more information on what each action computes
MOLINFOThis command is used to provide information on the molecules that are present in your system. More details STRUCTUREa file in pdb format containing a reference structure= ATOM 1 HH31 ACE 1 -9.105 -2.402 21.804 0.00 0.00
ATOM 2 CH3 ACE 1 -8.930 -3.352 22.308 0.00 0.00
ATOM 3 HH32 ACE 1 -9.504 -3.501 23.223 0.00 0.00
ATOM 4 HH33 ACE 1 -9.067 -4.173 21.604 0.00 0.00
ATOM 5 C ACE 1 -7.450 -3.303 22.659 0.00 0.00
...
ATOM 129 CH3 NME 14 -1.690 -7.089 20.487 0.00 0.00
ATOM 130 HH31 NME 14 -1.873 -6.629 21.459 0.00 0.00
ATOM 131 HH32 NME 14 -1.152 -6.284 19.986 0.00 0.00
ATOM 132 HH33 NME 14 -1.135 -8.020 20.602 0.00 0.00
END
The MOLINFO action with label calculates somethingWHOLEMOLECULESThis action is used to rebuild molecules that can become split by the periodic boundary conditions. More details ENTITY0the atoms that make up a molecule that you wish to align=1-100
alphaThe ALPHARMSD action with label alpha calculates the following quantities: Quantity | Type | Description |
alpha | vector | if LESS_THAN is present the RMSD distance between each residue and the ideal alpha helix. If LESS_THAN is not present the number of residue segments where the structure is similar to an alpha helix |
alpha_lessthan | scalar | the number blocks of residues that have an RMSD from the secondary structure that is less than the threshold |
: ALPHARMSDProbe the alpha helical content of a protein structure. This action is a shortcut. More details RESIDUESthis command is used to specify the set of residues that could conceivably form part of the secondary structure=all TYPE the manner in which RMSD alignment is performed=OPTIMAL LESS_THANcalculate the number of a residue segments that are within a certain target distance of this secondary structure type. Options for this keyword are explained in the documentation for LESS_THAN.={RATIONAL R_0=0.1 NN=8 MM=12}
# alpha: ALPHARMSD RESIDUES=all TYPE=OPTIMAL LESS_THAN={RATIONAL R_0=0.1 NN=8 MM=12}
alphaThe SECONDARY_STRUCTURE_RMSD action with label alpha calculates the following quantities: Quantity | Type | Description |
alpha | vector | a vector containing the rmsd distance between each of the residue segments and the reference structure |
: SECONDARY_STRUCTURE_RMSDCalclulate the distance between segments of a protein and a reference structure of interest More details BONDLENGTHthe length to use for bonds=0.17 ATOMSthis is the full list of atoms that we are investigating=7,9,11,15,16,17,19,21,25,26,27,29,31,35,36,37,39,41,45,46,47,49,51,55,56,57,59,61,65,66,67,69,71,75,76,77,79,81,85,86,87,89,91,95,96,97,99,101,105,106,107,109,111,115,116,117,119,121,125,126 TYPE the manner in which RMSD alignment is performed=OPTIMAL SEGMENT1this is the lists of atoms in the segment that are being considered=0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,17,18,19,20,21,22,23,24,25,26,27,28,29 SEGMENT2this is the lists of atoms in the segment that are being considered=5,6,7,8,9,10,11,12,13,14,15,16,17,18,19,20,21,22,23,24,25,26,27,28,29,30,31,32,33,34 SEGMENT3this is the lists of atoms in the segment that are being considered=10,11,12,13,14,15,16,17,18,19,20,21,22,23,24,25,26,27,28,29,30,31,32,33,34,35,36,37,38,39 SEGMENT4this is the lists of atoms in the segment that are being considered=15,16,17,18,19,20,21,22,23,24,25,26,27,28,29,30,31,32,33,34,35,36,37,38,39,40,41,42,43,44 SEGMENT5this is the lists of atoms in the segment that are being considered=20,21,22,23,24,25,26,27,28,29,30,31,32,33,34,35,36,37,38,39,40,41,42,43,44,45,46,47,48,49 STRUCTURE1the reference structure=0.733,0.519,5.298,1.763,0.81,4.301,3.166,0.543,4.881,1.527,-0.045,3.053,1.646,0.436,1.928,1.18,-1.312,3.254,0.924,-2.203,2.126,0.65,-3.626,2.626,-0.239,-1.711,1.261,-0.19,-1.815,0.032,-1.28,-1.172,1.891,-2.416,-0.661,1.127,-3.548,-0.217,2.056,-1.964,0.529,0.276,-2.364,0.659,-0.88,-1.13,1.391,0.856,-0.62,2.565,0.148,0.228,3.439,1.077,0.231,2.129,-1.032,0.179,2.733,-2.099,1.028,1.084,-0.833,1.872,0.593,-1.919,2.85,-0.462,-1.397,1.02,0.02,-3.049,1.317,0.227,-4.224,-0.051,-0.684,-2.696,-0.927,-1.261,-3.713,-1.933,-2.219,-3.074,-1.663,-0.171,-4.475,-1.916,-0.296,-5.673
alpha_ltThe LESS_THAN action with label alpha_lt calculates the following quantities: Quantity | Type | Description |
alpha_lt | vector | the vector obtained by doing an element-wise application of a function that is one if the input is less than a threshold to the input vectors |
: LESS_THANUse a switching function to determine how many of the input variables are less than a certain cutoff. More details ARGthe values input to this function=alpha SWITCHThis keyword is used if you want to employ an alternative to the continuous swiching function defined above={RATIONAL R_0=0.1 NN=8 MM=12}
alpha_lessthanThe SUM action with label alpha_lessthan calculates the following quantities: Quantity | Type | Description |
alpha_lessthan | scalar | the sum of all the elements in the input vector |
: SUMCalculate the sum of the arguments More details ARGthe values input to this function=alpha_lt PERIODICif the output of your function is periodic then you should specify the periodicity of the function=NO
# --- End of included input --- PRINTPrint quantities to a file. More details ARGthe labels of the values that you would like to print to the file=alpha.lessthan FILEthe name of the file on which to output these quantities=colvar
References
More information about how this action can be used is available in the following articles:
- F. Pietrucci, A. Laio, A Collective Variable for the Efficient Exploration of Protein Beta-Sheet Structures: Application to SH3 and GB1. Journal of Chemical Theory and Computation. 5, 2197–2201 (2009)
Syntax
The following table describes the keywords and options that can be used with this action
Keyword |
Type |
Default |
Description |
RESIDUES |
input |
none |
this command is used to specify the set of residues that could conceivably form part of the secondary structure |
TYPE |
compulsory |
DRMSD |
the manner in which RMSD alignment is performed |
D_0 |
compulsory |
0.0 |
The d_0 parameter of the switching function |
NN |
compulsory |
8 |
The n parameter of the switching function |
MM |
compulsory |
12 |
The m parameter of the switching function |
SERIAL |
optional |
false |
do the calculation in serial |
NOPBC |
optional |
false |
ignore the periodic boundary conditions |
VERBOSE |
optional |
false |
write a more detailed output |
LESS_THAN |
optional |
not used |
calculate the number of a residue segments that are within a certain target distance of this secondary structure type. Options for this keyword are explained in the documentation for LESS_THAN. |
R_0 |
optional |
not used |
The r_0 parameter of the switching function |