Tools

These tools may call external programs to perform each task. Please verify the installation of the necessary tool for each case.

Check the stereochemistry of protein residues

PDBTools.zetaFunction
zeta(r::Residue)

Computes the Cα chirality (zeta "virtual" torsion angle - Cα-N-C-Cβ).

Returns the torsion angle or NaN if the residue is not recognized a protein residue or if its a Gly residue. Expected values are 33.9 ± 3.5 degrees (for one standard deviation). Also see the zeta_check function.

Example

julia> using PDBTools

julia> protein = select(read_pdb(PDBTools.TESTPDB), "protein");

julia> residues = collect(eachresidue(protein));

julia> zeta(residues[1])
33.672016f0
source
PDBTools.zeta_checkFunction
zeta_check(r::Residue; nsigma=2)

Checks if the Cα chirality falls into expected ranges. See the zeta function for further information. The expected mean is 33.9 degrees with a standard deviation of 3.5 degrees. By default, nsigma=2, implying that the function returns true if the torsion falls within two standard deviations from the mean.

See: https://www.ebi.ac.uk/thornton-srv/software/PROCHECK/manual/manappa.html

Example

julia> using PDBTools

julia> protein = select(read_pdb(PDBTools.TESTPDB), "protein");

julia> residues = collect(eachresidue(protein));

julia> zeta_check(residues[1])
true
source

Add hydrogens with OpenBabel

PDBTools.add_hydrogens!Function
add_hydrogens!(atoms::AbstractVector{<:Atom}; pH=7.0, obabel="obabel", debug=false)

Add hydrogens to a PDB file using Open Babel.

Arguments

  • atoms::AbstractVector{<:Atom}: structure (usually PDB file of a protein) to add hydrogens to.
  • pH: the pH of the solution. Default is 7.0.
  • obabel: path to the obabel executable. Default is "obabel".
  • debug: if true, print the output message from obabel. Default is false.
Note

This function requires the installation of OpenBabel. Please cite the corresponding reference if using it.

Example

julia> using PDBTools

julia> atoms = read_pdb(PDBTools.TESTPDB, "protein and not element H");

julia> add_hydrogens!(atoms)
   Vector{Atom{Nothing}} with 1459 atoms with fields:
   index name resname chain   resnum  residue        x        y        z occup  beta model segname index_pdb
       1    N     ALA     A        1        1   -9.229  -14.861   -5.481  1.00  0.00     1       -         1
       2   CA     ALA     A        1        1   -8.483  -14.912   -6.726  1.00  0.00     1       -         2
       3   CB     ALA     A        1        1   -9.383  -14.465   -7.880  1.00  0.00     1       -         3
                                                       ⋮ 
    1457    H     THR     A      104      208    5.886  -10.722   -7.797  1.00  0.00     1       -      1457
    1458    H     THR     A      104      208    5.871  -10.612   -9.541  1.00  0.00     1       -      1458
    1459    H     THR     A      104      208    6.423  -12.076   -8.762  1.00  0.00     1       -      1459
source

Custom protein residue types

It is possible to add to the list of protein residues, custom residue types. This can be done by simply adding to the PDBTools.protein_residues dictionary of residues a new PDBTools.ProteinResidue entry. For example, here we create a new resiude type NEW with the same properties of an ALA residue. To remove all custom protein residues, use remove_custom_protein_residues!().

julia> using PDBTools

julia> remove_custom_protein_residues!();

julia> add_protein_residue!("NEW", PDBTools.protein_residues["ALA"])
PDBTools.ProteinResidue("NEW", "ALA", "A", "Aliphatic", false, true, 71.037114, 71.0779, 0, true)

julia> atom = Atom(resname="NEW");

julia> isprotein(atom)
true

julia> remove_custom_protein_residues!();

Here we repeteadly call remove_custom_residues!() to guarantee the proper execution of the test codes, without any custom residues in the list of protein residues.

PDBTools.add_protein_residue!Function
add_protein_residue!(resname::String, reference_residue::PDBTools.ProteinResidue)

Function to add a custom protein residue to the list of protein residues. The function will return the ProteinResidue object that was added. To remove all custom protein residues use remove_custom_protein_residues!().

Example

julia> using PDBTools

julia> remove_custom_protein_residues!();

julia> add_protein_residue!("sA", PDBTools.protein_residues["ALA"])
PDBTools.ProteinResidue("sA", "ALA", "A", "Aliphatic", false, true, 71.037114, 71.0779, 0, true)

julia> isprotein(Atom(resname="sA"))
true

julia> remove_custom_protein_residues!(); # clean up

Here we repeteadly call remove_custom_residues!() to guarantee the proper execution of the test codes, without any custom residues in the list of protein residues.

source
PDBTools.remove_custom_protein_residues!Function
remove_custom_protein_residues!()

Function to remove all custom protein residues from the list of protein residues.

Example

julia> using PDBTools

julia> remove_custom_protein_residues!(); # clean up

julia> add_protein_residue!("sA", PDBTools.protein_residues["ALA"])
PDBTools.ProteinResidue("sA", "ALA", "A", "Aliphatic", false, true, 71.037114, 71.0779, 0, true)

julia> isprotein(Atom(resname="sA"))
true

julia> remove_custom_protein_residues!();

julia> isprotein(Atom(resname="sA"))
false

Here we repeteadly call remove_custom_residues!() to guarantee the proper execution of the test codes, without any custom residues in the list of protein residues.

source

The SIRAH force-field residues and element types

Conveniencie functions can be created to add sets of new types of residues and atom types to the list of residues and elements. This is illustrated in the custom_types.jl file of the source code, in this case for the residues and atom types of the SIRAH force field for Coarse-Grained protein simulations.

With those definitions, adding all SIRAH protein residue types and element names can be done with:

julia> using PDBTools 

julia> remove_custom_protein_residues!(); remove_custom_elements!();

julia> custom_protein_residues!(SIRAH)
┌ Warning: 
│ 
│     Residue `sX` will be interpreted as bridged Cysteine.
│ 
└ @ PDBTools

julia> custom_elements!(SIRAH)
┌ Warning:
│
│     The element masses are not the coarse-grained ones. This must be fixed in the future.
│
└ @ PDBTools

julia> sirah_pdb = read_pdb(PDBTools.SIRAHPDB);

julia> resname.(eachresidue(sirah_pdb))
5-element Vector{InlineStrings.String7}:
 "sI"
 "sR"
 "sX"
 "sI"
 "sG"

julia> getseq(sirah_pdb)
5-element Vector{String}:
 "I"
 "R"
 "C"
 "I"
 "G"

julia> all(isprotein.(sirah_pdb))
true

julia> remove_custom_protein_residues!(); remove_custom_elements!();

Note that the residue names of the SIRAH force-field are non-standard (sI, sR, etc.), but the sequence is properly retrieved with standard one-letter codes, and all the atoms of the structure are recognized as being "protein" atoms.

Here we repeteadly call remove_custom_residues!() and remove_custom_elements!() to guarantee the proper execution of the test codes.

Move atoms and center of mass

The center_of_mass function can be used to compute the center of mass of set of atoms, and the moveto! function can be used to move the center of mass of the atoms to the origin (by default) or to a specified position:

MolSimToolkitShared.center_of_massFunction
center_of_mass(atoms::AbstractVector{<:Atom})

Calculate the center of mass of the atoms.

Example

julia> using PDBTools

julia> atoms = read_pdb(PDBTools.SMALLPDB);

julia> center_of_mass(atoms)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
  -5.584422772707132
 -13.110413081059928
  -7.139970851058855
source
PDBTools.moveto!Function
moveto!(atoms::AbstractVector{<:Atom}; center::AbstractVector{<:Real}=SVector(0.0, 0.0, 0.0))

Move the center of mass of the atoms to the specified center position, which defaults to the origin.

Example

julia> using PDBTools

julia> atoms = read_pdb(PDBTools.SMALLPDB);

julia> center_of_mass(atoms)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
  -5.584422772707132
 -13.110413081059928
  -7.139970851058855

julia> moveto!(atoms; center = [1.0, 2.0, 3.0]);

julia> center_of_mass(atoms)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
 1.0000000263619948
 1.9999999934852166
 2.9999999509668918
source