BioMakie API
BioMakie.UniProtData
BioMakie.atomcolors
BioMakie.atomradii
BioMakie.atomradius
BioMakie.atomsizes
BioMakie.backbonebonds
BioMakie.bestoccupants
BioMakie.bondshape
BioMakie.bondshapes
BioMakie.covalentbonds
BioMakie.distancebonds
BioMakie.firstlabel
BioMakie.getbonds
BioMakie.getbonds
BioMakie.getbonds
BioMakie.getinspectorlabel
BioMakie.getuniprotdata
BioMakie.msavalues
BioMakie.notwater
BioMakie.plotmsa
BioMakie.plotmsa!
BioMakie.plotstruc
BioMakie.plotstruc!
BioMakie.plottingdata
BioMakie.plottingdata
BioMakie.readuniprotdata
BioMakie.rescolors
BioMakie.showsummary
BioMakie.sidechainbonds
MakieCore.heatmap
MakieCore.heatmap
MakieCore.heatmap
MakieCore.heatmap!
MakieCore.heatmap!
MakieCore.heatmap!
BioMakie.UniProtData
— TypeUniProtData::DataType
A struct containing all the information from a UniProt entry.
General information
accession
id
protinfo
gene
gene_synonyms
secondary_accession
sequence
organism
datainfo
Features
molecule_processing
domains_and_sites
structural
ptm
sequence_information
mutagenesis
variants
topology
other_features
Comments
func
catalytic_activity
subunit
subcellular_location
interaction
tissue_specificity
polymorphism
allergen
web_resource
similarity
miscellaneous
other_comments
Database references (EMBL, PDB, etc.)
dbrefs
BioMakie.atomcolors
— Methodatomcolors( atoms )
Get a Vector of colors for the atoms. To see all default element and amino acid colorschemes, use getbiocolors()
. Keyword argument colors
takes a Dict which maps element to color. ("C" => :red)
This function uses 'MIToS.PDB.bestoccupancy' or 'defaultatom' to ensure only one position per atom.
Keyword Arguments:
- colors –- elecolors | Options - elecolors, aquacolors
BioMakie.atomradii
— Methodatomradii( atoms )
Collect atom radii based on element for plotting.
Keyword Arguments:
- radiustype –- :ballandstick | Options - :cov, :covalent, :vdw, :vanderwaals, :bas, :ballandstick, :spacefilling
BioMakie.atomradius
— Methodatomradius( atom )
Collect atom radius based on element for plotting.
Keyword Arguments:
- radiustype –- :ballandstick | Options - :cov, :covalent, :vdw, :vanderwaals, :bas, :ballandstick, :spacefilling
BioMakie.atomsizes
— Methodatomsizes( atms )
Get a Vector of sizes for the atoms from a BioStructures.StructuralElementOrList.
This function uses 'MIToS.PDB.bestoccupancy' or 'defaultatom' to ensure only one position per atom.
Keyword Arguments:
- radiustype –- :ballandstick | Options - :cov, :covalent, :vdw, :vanderwaals, :bas, :ballandstick, :spacefilling
BioMakie.backbonebonds
— Methodbackbonebonds( chn::BioStructures.Chain ) -> BitMatrix
Returns a matrix of backbone bonds in chn
, where Mat[i,j] = 1 if atoms i and j are bonded.
Keyword Arguments:
- cutoff –––––- 1.6 # distance cutoff for bonds
BioMakie.bestoccupants
— Methodbestoccupants( residues::Vector{MIToS.PDB.PDBResidue} )
Get an OrderedDict of the best occupancy atoms for each residue.
BioMakie.bondshape
— Methodbondshape( twoatoms )
bondshape( twopoints )
Returns a (mesh) cylinder between two atoms or atomic coordinates.
Keyword Arguments:
- bondwidth ––––––- 0.2
BioMakie.bondshapes
— Methodbondshapes( structure )
bondshapes( residues )
bondshapes( coordinates )
bondshapes( structure, bondmatrix )
bondshapes( residues, bondmatrix )
bondshapes( coordinates, bondmatrix )
Returns a (mesh) cylinder between two atoms or points.
Keyword Arguments:
- algo ––––––––– :knowledgebased | :distance, :covalent # unless bondmatrix is given
- distance ––––––– 1.9 # unless bondmatrix is given
- bondwidth ––––––- 0.2
BioMakie.covalentbonds
— Methodcovalentbonds( atms ) -> BitMatrix
Returns a matrix of all bonds in atms
, where Mat[i,j] = 1 if atoms i and j are bonded.
This function uses 'bestoccupancy' or 'defaultatom' to ensure only one position per atom.
Keyword Arguments:
- extradistance –– 0.14 # fudge factor for better inclusion
- H –––––––– true # include bonds with hydrogen atoms
- disulfides –––- false # include disulfide bonds
BioMakie.distancebonds
— Methoddistancebonds( atms ) -> BitMatrix
Returns a matrix of all bonds in atms
, where Mat[i,j] = 1 if atoms i and j are bonded.
This function uses 'bestoccupancy' or 'defaultatom' to ensure only one position per atom.
Keyword Arguments:
- cutoff –––––- 1.9 # distance cutoff for bonds between heavy atoms
- hydrogencutoff –- 1.14 # distance cutoff for bonds with hydrogen atoms
- H –––––––– true # include bonds with hydrogen atoms
- disulfides –––- false # include disulfide bonds
BioMakie.firstlabel
— Methodfirstlabel( inspectorfunc::Function )
Show an example of the inspector label function looks like. The position p
will not be available to this function, so it will be set to nothing
.
BioMakie.getbonds
— Methodgetbonds( chn::BioStructures.Chain, selectors... ) -> BitMatrix
getbonds( modl::BioStructures.Model, selectors... ) -> BitMatrix
getbonds( struc::BioStructures.MolecularStructure, selectors... ) -> BitMatrix
Returns a matrix of all bonds in chn
, where Mat[i,j] = 1 if atoms i and j are bonded.
This function uses 'bestoccupancy' or 'defaultatom' to ensure only one position per atom.
Keyword Arguments:
- algo ––––––- :knowledgebased # (:distance, :covalent) algorithm to find bonds
- H –––––––– true # include bonds with hydrogen atoms
- cutoff –––––- 1.9 # distance cutoff for bonds between heavy atoms
- extradistance –– 0.14 # fudge factor for better inclusion
- disulfides –––- false # include disulfide bonds
BioMakie.getbonds
— Methodgetbonds( residues ) -> BitMatrix
Returns a matrix of all bonds in residues::Vector{MIToS.PDB.PDBResidue}
, where Mat[i,j] = 1 if atoms i and j are bonded.
Keyword Arguments:
- algo ––––––- :knowledgebased # (:distance, :covalent) algorithm to find bonds
- H –––––––– true # include bonds with hydrogen atoms
- cutoff –––––- 1.9 # distance cutoff for bonds between heavy atoms
- extradistance –– 0.14 # fudge factor for better inclusion
- disulfides –––- false # include disulfide bonds
BioMakie.getbonds
— Methodgetbonds( coords ) -> BitMatrix
Returns a matrix of all bonds using a N x 3 coordinates matrix. Uses a plain cutoff distance with algo option :distance. This is not recommended as it can lead to incorrect results since different atoms have different bond lengths and radii.
Keyword Arguments:
- algo ––––––- :distance # algorithm to find bonds
- H –––––––– true # include bonds with hydrogen atoms
- cutoff –––––- 1.9 # distance cutoff for bonds between heavy atoms
- extradistance –– 0.14 # fudge factor for better inclusion
- disulfides –––- false # include disulfide bonds
BioMakie.getinspectorlabel
— Methodgetinspectorlabel( structure )
getinspectorlabel( residues )
getinspectorlabel( atom )
Get the inspector label function for plotting a 'StructuralElementOrList'.
This function uses 'MIToS.PDB.bestoccupancy' or 'defaultatom' to ensure only one position per atom.
BioMakie.getuniprotdata
— Functiongetuniprotdata(accession, filename=nothing; include_refs = false)
Downloads and reads a UniProt JSON file and returns a UniProtData struct.
Keyword Arguments:
- include_refs::Bool = false Whether to include all the database references (EMBL, PDB, etc.) in the struct. can be very large, so it is
false
by default.
BioMakie.msavalues
— Functionmsavalues( msa::AbstractMatrix, resdict::AbstractDict )::Matrix{Real}
Returns a matrix of numbers according to the given dictionary, where keys are residue letters and values are numbers. This matrix is used as input for plotmsa
for the heatmap colors of the residue positions.
Default values for residue letters are from Kidera Factor values. From: Kenta Nakai, Akinori Kidera, Minoru Kanehisa, Cluster analysis of amino acid indices for prediction of protein structure and function, Protein Engineering, Design and Selection, Volume 2, Issue 2, July 1988, Pages 93–100, https://doi.org/10.1093/protein/2.2.93
kf 2 is Kidera Factor 2 (size/volume-related). The KF dictionary is in utils.jl
, or you can look at the kideradict
variable.
Keyword Arguments:
- resdict –––- kideradict by default, alternatively give a Dict{String,Vector{Float}}
- kf –––––– 2 by default, alternatively give an integer from 1:10
BioMakie.notwater
— Methodnotwater( residues::Vector{MIToS.PDB.PDBResidue} )
Filter out water molecules from a Vector of residues.
BioMakie.plotmsa!
— Methodplotmsa!( fig, msa )
Plot a multiple sequence alignment (MSA) into a Figure.
Example
fig = Figure(resolution = (1100, 400))
plotmsa!( fig::Figure, msa::T; kwargs... ) where {T<:Union{MSA.AbstractMultipleSequenceAlignment,
Vector{Tuple{String,String}},
Vector{FASTX.FASTA.Record}}}
Keyword Arguments:
- sheetsize ––- [40,20]
- gridposition – (1,1)
- markersize –– 12
- colorscheme –- :buda
- markercolor –- :black
- xticklabelsize - 11
- yticklabelsize - 11
- resolution ––- (700,300)
- kwargs... # forwarded to scatter plot
BioMakie.plotmsa
— Methodplotmsa( msa )
plotmsa( plotdata )
Plot a multiple sequence alignment (MSA). Returns a Figure, or a Figure and Observables for interaction.
Examples
MIToS.Pfam.downloadpfam("PF00062") # download PF00062 MSA
msa = MIToS.MSA.read_file("PF00062.stockholm.gz", Stockholm,
generatemapping =true, useidcoordinates=true)
plotmsa( msa; kwargs... )
Keyword Arguments:
- figresolution ––- (1000,350) # because
resolution
applies to the MSA plot - sheetsize ––––- [40,20]
- gridposition ––– (1,1:3)
- colorscheme –––- :buda
- markersize –––– 12
- markercolor –––- :black
- xticklabelsize –– 11
- yticklabelsize –– 11
- kwargs... # forwarded to scatter plot
BioMakie.plotstruc!
— Methodplotstruc!( fig, structure )
plotstruc!( gridposition, structure )
plotstruc!( fig, plotdata )
plotstruc!( gridposition, plotdata )
Plot a protein structure(/chain/residues/atoms) into a Figure.
Examples
fig = Figure()
using MIToS.PDB
pdbfile = MIToS.PDB.downloadpdb("2vb1", format=PDBFile)
struc = MIToS.PDB.read_file(pdbfile, PDBFile) |> Observable
strucplot = plotstruc!(fig, struc)
chain_A = select_residues(struc, model="1", chain="A", group="ATOM")
strucplot = plotstruc!(fig, chain_A)
chnatms = select_atoms(struc, model="1", chain="A", group="ATOM")
strucplot = plotstruc!(fig, chnatms)
-------------------------
using BioStructures
struc = retrievepdb("2vb1", dir = "data/") |> Observable
strucplot = plotstruc!(fig, struc)
struc = read("data/2vb1_mutant1.pdb", BioStructures.PDBFormat) |> Observable
strucplot = plotstruc!(fig, struc)
chain_A = retrievepdb("2vb1", dir = "data/")["A"] |> Observable
strucplot = plotstruc!(fig, chain_A)
Keyword Arguments:
- resolution ––- (600,600)
- gridposition –- (1,1) # if an MSA is already plotted, (2,1:3) works well
- plottype –––- :ballandstick, :covalent, or :spacefilling
- atomcolors ––- elecolors, others in
getbiocolors()
, or provide a Dict like: "N" => :blue - markersize ––- 0.0
- markerscale –– 1.0
- bondtype –––- :knowledgebased, :covalent, or :distance
- distance –––- 1.9 # distance cutoff for covalent bonds
- inspectorlabel - :default, or define your own function like: (self, i, p) -> "atom: ... coords: ..."
- water ––––– false # show water molecules
- kwargs... ––– keyword arguments passed to the atom
meshscatter
BioMakie.plotstruc
— Methodplotstruc( structure )
plotstruc( residues )
plotstruc( plotdata )
Create and return a Makie Figure for a protein structural element.
Examples
using MIToS.PDB
pdbfile = MIToS.PDB.downloadpdb("2vb1", format=PDBFile)
struc = MIToS.PDB.read_file(pdbfile, PDBFile) |> Observable
strucplot = plotstruc(struc)
chain_A = select_residues(struc, model="1", chain="A", group="ATOM")
strucplot = plotstruc(chain_A)
chnatms = select_atoms(struc, model="1", chain="A", group="ATOM")
strucplot = plotstruc(chnatms)
-------------------------
using BioStructures
struc = retrievepdb("2vb1", dir = "data/") |> Observable
strucplot = plotstruc(struc)
struc = read("data/2vb1_mutant1.pdb", BioStructures.PDBFormat) |> Observable
strucplot = plotstruc(struc)
chain_A = retrievepdb("2hhb", dir = "data/")["A"] |> Observable
strucplot = plotstruc(chain_A)
Keyword Arguments:
- figresolution – (600,600) # because
resolution
applies to the plot - resolution ––- (600,600)
- gridposition –- (1,1) # if an MSA is already plotted, (2,1:3) works well
- plottype –––- :ballandstick, :covalent, or :spacefilling
- atomcolors ––- elecolors, others in
getbiocolors()
, or provide a Dict like: "N" => :blue - markersize ––- 0.0
- markerscale –– 1.0
- bondtype –––- :knowledgebased, :covalent, or :distance
- distance –––- 1.9 # distance cutoff for covalent bonds
- inspectorlabel - :default, or define your own function like: (self, i, p) -> "atom: ... coords: ..."
- water ––––– false # show water molecules
- kwargs... ––– keyword arguments passed to the atom
meshscatter
BioMakie.plottingdata
— Methodplottingdata( structure )
plottingdata( residues )
plottingdata( atoms )
This function returns an OrderedDict of the main data used for plotting. This function uses 'MIToS.PDB.bestoccupancy' or 'defaultatom' to ensure only one position per atom. By default the kwarg 'water' is set to false, so water molecules are not included.
Returns:
OrderedDict(:atoms => ...,
:coords => ...,
:colors => ...,
:sizes => ...,
:bonds => ...)
Keyword Arguments:
- colors –––- elecolors | Options - elecolors, aquacolors, shapelycolors, maecolors
- radiustype –- :ballandstick | Options - :cov, :covalent, :vdw, :vanderwaals, :bas, :ballandstick, :spacefilling
- water –––– false | Options - true, false
BioMakie.plottingdata
— Methodplottingdata( msa )
Collects data for plotting (residue string matrix, matrix heatmap values, x labels, and y labels) from a multiple sequence alignment (MSA) object.
The MSA object can be a:
AbstractMultipleSequenceAlignment
from MIToS.MSA,- vector of tuples 'Vector{Tuple{String,String}}' from FastaIO,
- vector of FASTA records 'Vector{FASTX.FASTA.Record}' from FASTX.
BioMakie.readuniprotdata
— Methodreaduniprotdata(jsonfile; include_refs = false)
Reads a UniProt JSON file and returns a UniProtData struct.
Keyword Arguments:
- include_refs::Bool = false Whether to include all the database references (EMBL, PDB, etc.) in the struct. can be very large, so it is
false
by default.
BioMakie.rescolors
— Methodrescolors( residues )
Get a Vector of colors for the atoms. To see all default element and amino acid colorschemes, use getbiocolors()
. Keyword argument colors
takes a Dict which maps residue to color. ("C" => :red)
This function uses 'MIToS.PDB.bestoccupancy' or 'defaultatom' to ensure only one position per atom.
Keyword Arguments:
- colors –- elecolors | Options - elecolors, aquacolors, shapelycolors, maecolors
BioMakie.showsummary
— Methodshowsummary(pdata)
Prints some of the most important information from a UniProtData object.
BioMakie.sidechainbonds
— Methodsidechainbonds( res::BioStructures.AbstractResidue, selectors... ) -> BitMatrix
Returns a matrix of sidechain bonds in res
, where Mat[i,j] = 1 if atoms i and j are bonded.
This function uses 'bestoccupancy' or 'defaultatom' to ensure only one position per atom.
Keyword Arguments:
- algo ––––––- :knowledgebased # (:distance, :covalent) algorithm to find bonds
- H –––––––– true # include bonds with hydrogen atoms
- cutoff –––––- 1.9 # distance cutoff for bonds between heavy atoms
- extradistance –– 0.14 # fudge factor for better inclusion
MakieCore.heatmap!
— Methodheatmap!( fig, dmap; kwargs... )
Plot a MIToS distance or contact map.
Example
fig = Figure()
using MIToS.PDB
pdbfile = MIToS.PDB.downloadpdb("1IVO", format=PDBFile)
residues_1ivo = read_file(pdbfile, PDBFile)
pdb = select_residues(residues_1ivo, model="1", chain="A", group="ATOM")
dmap = MIToS.PDB.distance(pdb, criteria="All") # distance map
cmap = contact(pdb, 8.0, criteria="CB") # contact map
heatmap!(fig, dmap)
Keyword Arguments:
- xlabel ––––––––- "Item 2"
- ylabel ––––––––- "Item 1"
- colormap –––––––- :viridis
- kwargs... ––––––– additional keyword arguments to pass to heatmap
MakieCore.heatmap!
— Methodheatmap!( fig, cmap; kwargs... )
Plot a BioStructures contact map.
Example
fig = Figure()
using BioStructures
struc = retrievepdb("1IVO")[1]
cbetas_A = collectatoms(struc["A"], cbetaselector)
cbetas_B = collectatoms(struc["B"], cbetaselector)
cmap = ContactMap(cbetas_A, cbetas_B)
heatmap!(fig, cmap)
Keyword Arguments:
- xlabel ––––––––- "Item 2"
- ylabel ––––––––- "Item 1"
- colormap –––––––- :ice
- kwargs... ––––––– Keyword arguments to pass to heatmap
MakieCore.heatmap!
— Methodheatmap( dmap; kwargs... )
Plot a BioStructures distance map.
Example
using BioStructures
struc = retrievepdb("1IVO")[1]
cbetas_A = collectatoms(struc["A"], cbetaselector)
cbetas_B = collectatoms(struc["B"], cbetaselector)
dmap = DistanceMap(cbetas_A, cbetas_B)
heatmap(dmap)
Keyword Arguments:
- xlabel ––––––––- "Item 2"
- ylabel ––––––––- "Item 1"
- colormap –––––––- :viridis
- kwargs... ––––––– additional keyword arguments to pass to heatmap
MakieCore.heatmap
— Methodheatmap( dmap; kwargs... )
Plot a MIToS distance or contact map.
Example
using MIToS.PDB
pdbfile = MIToS.PDB.downloadpdb("1IVO", format=PDBFile)
residues_1ivo = read_file(pdbfile, PDBFile)
pdb = select_residues(residues_1ivo, model="1", chain="A", group="ATOM")
dmap = MIToS.PDB.distance(pdb, criteria="All") # distance map
cmap = contact(pdb, 8.0, criteria="CB") # contact map
heatmap(dmap)
Keyword Arguments:
- xlabel ––––––––- "Item 2"
- ylabel ––––––––- "Item 1"
- colormap –––––––- :viridis
- kwargs... ––––––– additional keyword arguments to pass to heatmap
MakieCore.heatmap
— Methodheatmap( cmap; kwargs... )
Plot a BioStructures contact map.
Example
using BioStructures
struc = retrievepdb("1IVO")[1]
cbetas_A = collectatoms(struc["A"], cbetaselector)
cbetas_B = collectatoms(struc["B"], cbetaselector)
cmap = ContactMap(cbetas_A, cbetas_B)
heatmap(cmap)
Keyword Arguments:
- xlabel ––––––––- "Item 2"
- ylabel ––––––––- "Item 1"
- colormap –––––––- :ice
- kwargs... ––––––– additional keyword arguments to pass to heatmap
MakieCore.heatmap
— Methodheatmap( dmap; kwargs... )
Plot a BioStructures distance map.
Example
using BioStructures
struc = retrievepdb("1IVO")[1]
cbetas_A = collectatoms(struc["A"], cbetaselector)
cbetas_B = collectatoms(struc["B"], cbetaselector)
dmap = DistanceMap(cbetas_A, cbetas_B)
heatmap(dmap)
Keyword Arguments:
- xlabel ––––––––- "Item 2"
- ylabel ––––––––- "Item 1"
- colormap –––––––- :viridis
- kwargs... ––––––– additional keyword arguments to pass to heatmap