Help entries

These entries can be accessed from the Julia REPL by typing ?, for example,

julia> ? mass
search: mass mapslices MathConstants makedocs set_zero_subnormals get_zero_subnormals mutable struct

  mass(name::String or atom::Atom or Vector{Atom})

  Returns the mass of an atom given its name, or Atom structure, or the total mass of a vector of Atoms.

  Example
  –––––––––

  julia> atoms = [ Atom(name="NT3"), Atom(name="CA") ];
  
  julia> mass(atoms[1])
  14.0067
  
  julia> mass("CA")
  12.011
  
  julia> mass(atoms)
  26.017699999999998
PDBTools.AtomType
Atom::DataType

Structure that contains the atom properties. It is mutable, so it can be edited. Fields:

mutable struct Atom
    index::Int # The sequential index of the atoms in the file
    index_pdb::Int # The index as written in the PDB file (might be anything)
    name::String # Atom name
    resname::String # Residue name
    chain::String # Chain identifier
    resnum::Int # Number of residue as written in PDB file
    residue::Int # Sequential residue (molecule) number in file
    x::Float64 # x coordinate
    y::Float64 # y coordinate
    z::Float64 # z coordinate
    beta::Float64 # temperature factor
    occup::Float64 # occupancy
    model::Int # model number
    segname::String # Segment name (cols 73:76)
    pdb_element::String # Element symbol string (cols 77:78)
    charge::String # Charge (cols: 79:80)
    custom::Dict{Symbol, Any} # Custom fields
end

Example

julia> using PDBTools

julia> atoms = readPDB(PDBTools.SMALLPDB)
   Array{Atoms,1} with 35 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  0.00  0.00     1    PROT         1
       2 1HT1     ALA     A        1        1  -10.048  -15.427   -5.569  0.00  0.00     1    PROT         2
       3  HT2     ALA     A        1        1   -9.488  -13.913   -5.295  0.00  0.00     1    PROT         3
                                                       ⋮
      33  OD2     ASP     A        3        3   -6.974  -11.289   -9.300  1.00  0.00     1    PROT        33
      34    C     ASP     A        3        3   -2.626  -10.480   -7.749  1.00  0.00     1    PROT        34
      35    O     ASP     A        3        3   -1.940  -10.014   -8.658  1.00  0.00     1    PROT        35

julia> resname(atoms[1])
"ALA"

julia> chain(atoms[1])
"A"

julia> element(atoms[1])
"N"

julia> mass(atoms[1])
14.0067

julia> position(atoms[1])
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
  -9.229
 -14.861
  -5.481

The pdb_element and charge fields, which are frequently left empty in PDB files, are not printed. The direct access to the fields is considered part of the interface.

Custom fields can be set on Atom construction with the custom keyword argument, which receives a Dict{Symbol,Any} as parameter. They can be retrieved with the custom_field function or, if the custom field names does not overlap with an existing field, with the dot syntax. Requires PDBTools > 0.14.3.

Example

julia> using PDBTools

julia> atom = Atom(index = 0; custom=Dict(:c => "c", :index => 1));

julia> atom.c
"c"

julia> index(atom)
0

julia> custom_field(atom, :index)
1
source
PDBTools.ResidueType
Residue(atoms::AbstractVector{Atom}, range::UnitRange{Int})

Residue data structure. It contains two fields: atoms which is a vector of Atom elements, and range, which indicates which atoms of the atoms vector compose the residue.

The Residue structure carries the properties of the residue or molecule of the atoms it contains, but it does not copy the original vector of atoms, only the residue meta data for each residue.

Example

julia> pdb = wget("1LBD");

julia> residues = collect(eachresidue(pdb))
   Array{Residue,1} with 238 residues.

julia> resnum.(residues[1:3])
3-element Vector{Int64}:
 225
 226
 227

julia> residues[5].chain
"A"

julia> residues[8].range
52:58
source
PDBTools.SelectType
Select

This structure acts a function when used within typical julia filtering functions, by converting a string selection into a call to query call.

Example

julia> using PDBTools

julia> atoms = readPDB(PDBTools.TESTPDB, "protein");

julia> findfirst(Select("name CA"), atoms)
5

julia> filter(Select("name CA and residue 1"), atoms)
   Array{Atoms,1} with 1 atoms with fields:
   index name resname chain   resnum  residue        x        y        z occup  beta model segname index_pdb
       5   CA     ALA     A        1        1   -8.483  -14.912   -6.726  1.00  0.00     1    PROT         5
source
PDBTools.SequenceType
Sequence

Wrapper for strings, or vectors of chars, strings, or residue names, to dispatch on functions that operate on amino acid sequences.

Example

julia> seq = ["Alanine", "Glutamic acid", "Glycine"];

julia> mass(Sequence(seq))
257.2432

julia> seq = "AEG";

julia> mass(Sequence(seq))
257.2432
source
AtomsBase.atomic_numberMethod
atomic_number(atom::Atom)

Returns the atomic number of an atom from its Atom structure.

Example

julia> using PDBTools

julia> at = Atom(name="NT3");

julia> atomic_number(at)
7
source
PDBTools._closestMethod
closest(x,y)

Computes the minimum distance between two sets of atoms and returns the indices of the atoms and their distance. Both vector of atoms or vectors of coordinates can be used as input.

Examples

julia> model = wget("1BSX");

julia> protein = select(model,"protein");

julia> ligand = select(model,"resname T3");

julia> closest(ligand,protein)
(43, 3684, 2.7775834820937417)

julia> ligand[43]
    4037   O1      T3     B        2      512  -22.568   81.625    3.159 36.59  1.00     1       -      4041

julia> closest(ligand[43],protein)
(1, 3684, 2.7775834820937417)

julia> x = coor(protein)
3994-element Vector{SVector{3, Float64}}:
 [52.884, 24.022, 35.587]
 [52.916, 24.598, 36.993]
 ⋮
 [-46.887, 86.925, 13.235]
 [-47.164, 83.593, 15.25]

julia> closest(ligand,x)
(43, 3684, 2.7775834820937417)
source
PDBTools.add_element!Method
add_element!(symbol::String, reference_element::PDBTools.Element; elements=PDBTools.elements)

Add a new element to the elements dictionary. If the element already exists, overwrite it.

To remove all custom elements, use remove_custom_elements!().

Example

julia> using PDBTools

julia> remove_custom_elements!(); # if any

julia> atoms = [ Atom(name="A1"), Atom(name="A2") ];

julia> add_element!("A1", PDBTools.elements["C"])
PDBTools.Element(:C, "C", "Carbon", 6, 12.011, true)

julia> add_element!("A2", PDBTools.elements["N"])
PDBTools.Element(:N, "N", "Nitrogen", 7, 14.0067, true)

julia> element(atoms[1])
"C"

julia> element(atoms[2])
"N"

julia> mass(atoms)
26.017699999999998

julia> remove_custom_elements!(); 

Here we repeteadly call remove_custom_elements!() to guarantee the proper execution of the test codes, without any custom elements predefined.

source
PDBTools.add_hydrogens!Method
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 = readPDB(PDBTools.TESTPDB, "protein and not element H");

julia> add_hydrogens!(atoms)
   Array{Atoms,1} 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
PDBTools.add_protein_residue!Method
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, false, 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.center_of_massMethod
center_of_mass(atoms::AbstractVector{<:Atom})

Calculate the center of mass of the atoms.

Example

julia> using PDBTools

julia> atoms = readPDB(PDBTools.SMALLPDB);

julia> center_of_mass(atoms)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
  -5.584422752942997
 -13.110413157869903
  -7.139970815730879
source
PDBTools.coorMethod
coor(atoms; selection)

Returns the coordinates of the atoms. The input may be one atom (type Atom), a vector of atoms, or a Residue. The coordinates are returned as a vector of static vectors (from StaticArrays), more specifically as a Vector{SVector{3,Float64}}.

Examples

julia> using PDBTools, StaticArrays 

julia> protein = wget("1LBD");

julia> coor(protein[1])
3-element SVector{3, Float64} with indices SOneTo(3):
 45.228
 84.358
 70.638

julia> coor(protein[1],as=SVector{3,Float32})
3-element SVector{3, Float32} with indices SOneTo(3):
 45.228
 84.358
 70.638

julia> coor(protein, "index <= 2")
2-element Vector{SVector{3, Float64}}:
 [45.228, 84.358, 70.638]
 [46.08, 83.165, 70.327]

julia> coor(protein, only = at -> at.resname == "ALA")
110-element Vector{SVector{3, Float64}}:
 [43.94, 81.982, 70.474]
 [43.02, 80.825, 70.455]
 [41.996, 80.878, 69.34]
 ⋮
 [-17.866, 84.088, 51.741]
 [-18.496, 83.942, 52.777]
 [-15.888, 82.583, 51.706]
  
julia> residues = collect(eachresidue(protein));

julia> coor(residues[1])
6-element Vector{SVector{3, Float64}}:
 [45.228, 84.358, 70.638]
 [46.08, 83.165, 70.327]
 [45.257, 81.872, 70.236]
 [45.823, 80.796, 69.974]
 [47.147, 82.98, 71.413]
 [46.541, 82.639, 72.662]
source
PDBTools.distanceMethod
distance(x,y)

Computes the minimum distance between two sets of atoms, between an atom and a set of atoms, or simply the distance between two atoms. The input may be a vector of Atoms, or the coordinates that are output of the coor function.

Examples

julia> model = wget("1BSX");

julia> protein = select(model,"protein");

julia> ligand = select(model,"resname T3");

julia> distance(protein,ligand)
2.7775834820937417

julia> distance(protein[1],ligand[3])
36.453551075306784

julia> distance(coor(ligand),protein)
2.7775834820937417
source
PDBTools.eachresidueMethod
eachresidue(atoms::AbstractVector{Atom})

Iterator for the residues (or molecules) of a selection.

Example

julia> atoms = wget("1LBD");

julia> length(eachresidue(atoms))
238

julia> for res in eachresidue(atoms)
         println(res)
       end
 Residue of name SER with 6 atoms.
   index name resname chain   resnum  residue        x        y        z  beta occup model segname index_pdb
       1    N     SER     A      225        1   45.228   84.358   70.638 67.05  1.00     1       -         1
       2   CA     SER     A      225        1   46.080   83.165   70.327 68.73  1.00     1       -         2
       3    C     SER     A      225        1   45.257   81.872   70.236 67.90  1.00     1       -         3
       4    O     SER     A      225        1   45.823   80.796   69.974 64.85  1.00     1       -         4
       5   CB     SER     A      225        1   47.147   82.980   71.413 70.79  1.00     1       -         5
       6   OG     SER     A      225        1   46.541   82.639   72.662 73.55  1.00     1       -         6

 Residue of name ALA with 5 atoms.
   index name resname chain   resnum  residue        x        y        z  beta occup model segname index_pdb
       7    N     ALA     A      226        2   43.940   81.982   70.474 67.09  1.00     1       -         7
       8   CA     ALA     A      226        2   43.020   80.825   70.455 63.69  1.00     1       -         8
       9    C     ALA     A      226        2   41.996   80.878   69.340 59.69  1.00     1       -         9
                                                      ...
source
PDBTools.edit!Method
edit!(atoms::Vector{Atom})

Opens a temporary PDB file in which the fields of the vector of atoms can be edited.

source
PDBTools.elementMethod
element(atom::Atom)

Returns the element symbol, as a string, of an atom given the Atom structure. If the pdb_element is empty or "X", the element is inferred from the atom name. Othwerwise, the pdb_element is returned.

Example

julia> using PDBTools

julia> at = Atom(name="NT3");

julia> element(at)
"N"
source
PDBTools.element_nameMethod
element_name(atom::Atom)

Returns the element name of an atom given its name, or Atom structure.

Example

julia> using PDBTools

julia> at = Atom(name="NT3");

julia> element_name(at)
"Nitrogen"
source
PDBTools.element_symbolMethod
element_symbol(atom::Atom)

Returns a symbol for element name of an atom given its name, or Atom structure.

Example

julia> using PDBTools 

julia> at = Atom(name="NT3");

julia> element_symbol(at)
:N
source
PDBTools.element_symbol_stringMethod
element_symbol_string(atom::Atom)

Returns a string with the symbol of the element, given the Atom structure.

Example

julia> using PDBTools 

julia> at = Atom(name="NT3");

julia> element_symbol_string(at)
"N"
source
PDBTools.formulaMethod
formula(atoms::AbstractVector{Atom})

Returns the molecular formula of the current selection.

Example

julia> using PDBTools

julia> pdb  = readPDB(PDBTools.TESTPDB, "residue 1"); # testing PDB file

julia> resname(pdb[1])
"ALA"

julia> formula(pdb)
H₇C₃N₁O₁
source
PDBTools.getseqMethod
getseq(Vector{Atom} or filename; selection, code)

Returns the sequence of aminoacids from the vector of atoms or file name. Selections may be applied. Code defines if the output will be a one-letter, three-letter or full-residue name array.

Example

julia> protein = wget("1LBD");

julia> getseq(protein,"residue < 3")
2-element Vector{String}:
 "S"
 "A"

julia> getseq(protein,"residue < 3", code=2)
2-element Vector{String}:
 "SER"
 "ALA"

julia> getseq(protein,"residue < 3",code=3)
2-element Vector{String}:
 "Serine"
 "Alanine"
source
PDBTools.massMethod
mass(atom::Atom)
mass(atoms::AbstractVector{<:Atoms})

Returns the mass of an atom given its name, or Atom structure, or the total mass of a vector of Atoms.

If a mass is defined as a custom field in the the Atom structure, it is returned. Otherwise, the mass is retrieved from the element mass as inferred from the atom name.

Example

julia> using PDBTools

julia> atoms = [ Atom(name="NT3"), Atom(name="CA") ];

julia> mass(atoms[1])
14.0067

julia> mass(atoms)
26.017699999999998
source
PDBTools.massMethod
mass(s::Sequence)

Returns the mass of a sequence of amino acids, given a Sequence struct type.

Examples

julia> seq = ["Alanine", "Glutamic acid", "Glycine"];

julia> mass(Sequence(seq))
257.2432

julia> seq = "AEG";

julia> mass(Sequence(seq))
257.2432

julia> seq = ["ALA", "GLU", "GLY"];

julia> mass(Sequence(seq))
257.2432
source
PDBTools.maxminMethod
maxmin(atoms::Vector{Atom}; selection)

Returns the maximum and minimum coordinates of an atom vector, and the length (maximum minus minimum) in each direction.

Example

julia> protein = wget("1LBD");

julia> maxmin(protein)
 
 Minimum atom coordinates: xmin = [-29.301, 57.178, 45.668]
 Maximum atom coordinates: xmax = [47.147, 99.383, 86.886]
 Length in each direction: xlength = [76.448, 42.205, 41.217999999999996]
source
PDBTools.moveto!Method
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 = readPDB(PDBTools.SMALLPDB);

julia> center_of_mass(atoms)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
  -5.584422752942997
 -13.110413157869903
  -7.139970815730879

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.0
 2.0000000000000036
 3.0000000000000018
source
PDBTools.oneletterMethod
oneletter(residue::Union{String,Char})

Function to return a one-letter residue code from the three letter code or residue name. The function is case-insensitive.

Examples

julia> oneletter("ALA")
"A"

julia> oneletter("Glutamic acid")
"E"
source
PDBTools.printatomMethod
printatom(atom::Atom)

Prints an Atom structure in a human-readable format, with a title line.

Example

julia> using PDBTools

julia> atoms = readPDB(PDBTools.TESTPDB, "protein and residue 2");

julia> printatom(atoms[1])
   index name resname chain   resnum  residue        x        y        z occup  beta model segname index_pdb
      13    N     CYS     A        2        2   -6.351  -14.461   -5.695  1.00  0.00     1    PROT        13

julia> atoms[1] # default show method
      13    N     CYS     A        2        2   -6.351  -14.461   -5.695  1.00  0.00     1    PROT        13
source
PDBTools.readPDBFunction
readPDB(pdbfile::String, selection::String)
readPDB(pdbfile::String; only::Function = all)

readPDB(pdbdata::IOBuffer, selection::String)
readPDB(pdbdata::IOBuffer; only::Function = all)

Reads a PDB file and stores the data in a vector of type Atom.

If a selection is provided, only the atoms matching the selection will be read. For example, resname ALA will select all the atoms in the residue ALA.

If the only function keyword is provided, only the atoms for which only(atom) is true will be read.

Examples

julia> protein = readPDB("../test/structure.pdb")
   Array{Atoms,1} with 62026 atoms with fields:
   index name resname chain   resnum  residue        x        y        z  beta occup model segname index_pdb
       1    N     ALA     A        1        1   -9.229  -14.861   -5.481  0.00  1.00     1    PROT         1
       2  HT1     ALA     A        1        1  -10.048  -15.427   -5.569  0.00  0.00     1    PROT         2
                                                       ⋮ 
   62025   H1    TIP3     C     9339    19638   13.218   -3.647  -34.453  0.00  1.00     1    WAT2     62025
   62026   H2    TIP3     C     9339    19638   12.618   -4.977  -34.303  0.00  1.00     1    WAT2     62026

julia> ALA = readPDB("../test/structure.pdb","resname ALA")
   Array{Atoms,1} with 72 atoms with fields:
   index name resname chain   resnum  residue        x        y        z  beta occup model segname index_pdb
       1    N     ALA     A        1        1   -9.229  -14.861   -5.481  0.00  1.00     1    PROT         1
       2  HT1     ALA     A        1        1  -10.048  -15.427   -5.569  0.00  0.00     1    PROT         2
                                                       ⋮ 
    1339    C     ALA     A       95       95   14.815   -3.057   -5.633  0.00  1.00     1    PROT      1339
    1340    O     ALA     A       95       95   14.862   -2.204   -6.518  0.00  1.00     1    PROT      1340

julia> ALA = readPDB("../test/structure.pdb", only = atom -> atom.resname == "ALA")
   Array{Atoms,1} with 72 atoms with fields:
   index name resname chain   resnum  residue        x        y        z  beta occup model segname index_pdb
       1    N     ALA     A        1        1   -9.229  -14.861   -5.481  0.00  1.00     1    PROT         1
       2  HT1     ALA     A        1        1  -10.048  -15.427   -5.569  0.00  0.00     1    PROT         2
                                                       ⋮ 
    1339    C     ALA     A       95       95   14.815   -3.057   -5.633  0.00  1.00     1    PROT      1339
    1340    O     ALA     A       95       95   14.862   -2.204   -6.518  0.00  1.00     1    PROT      1340
source
PDBTools.remove_custom_elements!Function
remove_custom_elements!()

Remove all custom elements from the elements dictionary.

Example

julia> using PDBTools

julia> remove_custom_elements!();

julia> add_element!("GN", PDBTools.elements["N"])
PDBTools.Element(:N, "N", "Nitrogen", 7, 14.0067, true)

julia> element(Atom(name="GN"))
"N"

julia> remove_custom_elements!();

julia> element(Atom(name="GN")) # returns `nothing`

Here we repeteadly call remove_custom_elements!() to guarantee the proper execution of the test codes, without any custom elements predefined.

source
PDBTools.remove_custom_protein_residues!Method
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, false, 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
PDBTools.residue_ticksMethod
residue_ticks(atoms::AbstractVector{<:Atom}; first=nothing, last=nothing, stride=1, oneletter=true)

Returns a tuple with residue numbers and residue names for the given atoms, to be used as tick labels in plots.

first and last optional keyword parameters are integers that refer to the residue numbers to be included. The stride option can be used to skip residues and declutter the tick labels.

If oneletter is false, three-letter residue codes are returned. Residues with unknown names will be named X or XXX.

Examples

julia> using PDBTools

julia> atoms = wget("1UBQ", "protein");

julia> residue_ticks(atoms; stride=10)
([1, 11, 21, 31, 41, 51, 61, 71], ["M1", "K11", "D21", "Q31", "Q41", "E51", "I61", "L71"])

The resulting tuple of residue numbers and labels can be used as xticks in Plots.plot, for example.

source
PDBTools.residuenameMethod
residuename(residue::Union{String,Char})

Function to return the long residue name from other residue codes. The function is case-insensitive.

Examples

julia> residuename("A")
"Alanine"

julia> residuename("Glu")
"Glutamic Acid"
source
PDBTools.resnameMethod
resname(residue::Union{String,Char})

Returns the residue name, given the one-letter code or residue name. Differently from threeletter, this function will return the force-field name if available in the list of protein residues.

Examples

julia> resname("ALA")
"ALA"

julia> resname("GLUP")
"GLUP"
source
PDBTools.select_with_vmdMethod
select_with_vmd(inputfile::String, selection::String; vmd="vmd", srcload=nothing)
select_with_vmd(atoms::AbstractVector{<:Atom}, selection::String; vmd="vmd", srcload=nothing)

Select atoms using vmd selection syntax, with vmd in background. The input can be a file or a list of atoms.

Returns a tuple with list of index (one-based) and atom names of the selection.

Function to return the selection from a input file (topology, coordinates, etc), by calling VMD in the background.

The srcload argument can be used to load a list of scripts before loading the input file, for example with macros to define custom selection keywords.

source
PDBTools.stoichiometryMethod
stoichiometry(atoms::AbstractVector{Atom})

Returns the stoichiometry of atom selection in a Formula structure.

Example

julia> using PDBTools

julia> pdb  = readPDB(PDBTools.TESTPDB, "water"); # testing PDB file

julia> stoichiometry(pdb)
H₂O₁
source
PDBTools.threeletterMethod
threeletter(residue::String)

Function to return the three-letter natural-amino acid residue code from the one-letter code or residue name. The function is case-insensitive.

Examples

julia> threeletter("A")
"ALA"

julia> threeletter("Aspartic acid")
"ASP"

julia> threeletter("HSD")
"HIS"
source
PDBTools.wgetMethod
wget(PDBid; selection)

Retrieves a PDB file from the protein data bank. Selections may be applied.

Example

julia> protein = wget("1LBD","chain A")
   Array{Atoms,1} with 1870 atoms with fields:
   index name resname chain   resnum  residue        x        y        z  beta occup model segname index_pdb
       1    N     SER     A      225        1   45.228   84.358   70.638 67.05  1.00     1       -         1
       2   CA     SER     A      225        1   46.080   83.165   70.327 68.73  1.00     1       -         2
       3    C     SER     A      225        1   45.257   81.872   70.236 67.90  1.00     1       -         3
                                                       ⋮ 
    1868  OG1     THR     A      462      238  -27.462   74.325   48.885 79.98  1.00     1       -      1868
    1869  CG2     THR     A      462      238  -27.063   71.965   49.222 78.62  1.00     1       -      1869
    1870  OXT     THR     A      462      238  -25.379   71.816   51.613 84.35  1.00     1       -      1870
source
PDBTools.writePDBMethod
writePDB(atoms::Vector{Atom}, filename, selection; header=:auto, footer=:auto)

Write a PDB file with the atoms in atoms to filename. The selection argument is a string that can be used to select a subset of the atoms in atoms. For example, writePDB(atoms, "test.pdb", "name CA").

The header and footer arguments can be used to add a header and footer to the PDB file. If header is :auto, then a header will be added with the number of atoms in atoms. If footer is :auto, then a footer will be added with the "END" keyword. Either can be set to nothing if no header or footer is desired.

source