MASHUOA / HiTMaP

An R package of High-resolution Informatics Toolbox for Maldi-imaging Proteomics
GNU General Public License v3.0
15 stars 12 forks source link
docker-image hitmap maldi-imaging-proteomics

HiT-MaP

– An R package of High-resolution Informatics Toolbox for Maldi-imaging Proteomics

Find our published research article on Nature Communications:

https://doi.org/10.1038/s41467-021-23461-w

https://doi.org/10.1038/s41467-021-23461-w
zenodo

Maintainer: George Guo george.guo@auckland.ac.nz

About us:

Mass Spectrometry Hub | University of Auckland

Cancer research theme | Garvan Institute of Medical Research

MSRC Schey lab | Vanderbilt University

Package installation

This is a tutorial for the use of HiTMaP (An R package of High-resolution Informatics Toolbox for Maldi-imaging Proteomics). User’s may run HiTMaP using Docker, or through R console, however Docker is recommended to avoid issues with package dependency.

Installation of docker image

HiTMaP has been encapsulated into a docker image. After a proper installation and configuration of Docker engine (Docker documentation), user’s can download the latest version of HiT-MaP by using the bash code as below.

docker pull mashuoa/hitmap

Tags of available docker images:

  1. mashuoa/hitmap:latest contains the stable build release (built from the Dockerfile at MASHUOA/hitmap_docker with the effort from John Reeves j.reeves@garvan.org.au).

  2. mashuoa/hitmap:natcomms contains the original version when this project been accepted (minor changes applied to enhance the multi-files cluster image rendering).

  3. mashuoa/hitmap:gui_latest contains the developing graphical user interface of HiTMaP. Please map the 3838 port to the container and access the GUI via http://localhost:3838/. We are happy to hear your voice regarding the High-RES IMS pre-processing, segmentation and annotation as well as their corresponding GUI configurations.

  4. We are able to supply a singularity template to the users who want to deploy the HiTMaP on an HPC server. This scripts also are available at the MASHUOA/hitmap/dockerfiles.

Setting up and running the docker container:

# For windows user's, run the image with a local user\Documents\expdata folder mapped to the docker container:
docker run --name hitmap -v %userprofile%\Documents\expdata:/root/expdata -a stdin -a stdout -i -t mashuoa/hitmap /bin/bash 
# For linux or mac user's, run the image with a local user/expdata folder mapped to the docker container:
docker run --name hitmap -v ~/expdata:/root/expdata -a stdin -a stdout -i -t mashuoa/hitmap /bin/bash 

#Run the R console
R

Revoke Docker terminal:

#use ctrl+d to exit the docker container shell 

#Restart the container and connect to the shell 
docker restart hitmap
docker container exec -it hitmap /bin/bash

Stop/remove docker container (warning: if no local disk is mapped to “~/expdata”, please backup your existing result files from the container before you remove it):

docker stop hitmap
docker rm hitmap

If you are using docker GUI, pull the docker image using the codes above and follow the image as below to setup the container. if you are using mashuoa/hitmap:shiny_server, please also map local host:3838 to the container (Ports -> local hosts -> 3838).

Docker GUI setting

Installation code for R console installation

The code below is used for an experienced R user to build a local R/HiTMaP running environment. Major dependencies to note:

#install the git package
install.packages("remotes")
install.packages("devtools")
install.packages("BiocManager")
#library(devtools)
library(remotes)
Sys.setenv("R_REMOTES_NO_ERRORS_FROM_WARNINGS" = "true")
options(install.packages.check.source = "no")
BiocManager::install(c( "XVector", "Biostrings", "KEGGREST","cleaver"),INSTALL_opts="-Wno-error")
BiocManager::install(c("EBImage","Rdisop"))
remotes::install_github("MASHUOA/HiTMaP",force=T,build_opts = c("--no-resave-data", "--no-manual","-Wno-error", "--no-build-vignettes"),configure.vars="CFLAGS= -O3 -Wall -mtune=native -march=native -Wno-error",ask = F)

#Update all dependencies
BiocManager::install(ask = F)
library(HiTMaP)

Codes for Linux OS building enviornment

Run the codes as below to enable the required components in Linux console.

function apt_install() {
    if ! dpkg -s "$@" >/dev/null 2>&1; then
        if [ "$(find /var/lib/apt/lists/* | wc -l)" = "0" ]; then
            apt-get update
        fi
        apt-get install -y --allow-downgrades --no-install-recommends "$@"
    fi
}

apt_install \
    sudo \
    gdebi-core \
    libcairo2=1.18.0-1+b1 \
    libcairo-script-interpreter2=1.18.0-1+b1 \
    lsb-release \
    libcurl4-openssl-dev \
    libcairo2-dev \
    libxt-dev \
    xtail \
    wget \
    default-jdk \
    libxml2-dev \
    libssl-dev \
    libudunits2-dev \
    librsvg2-dev \
    libmagick++-dev \
    r-cran-ncdf4 \
    libz-dev \
    libnss-winbind \
    winbind \
    dirmngr \
    gnupg \
    apt-transport-https \
    ca-certificates \
    software-properties-common \
    libfftw3-dev \
    texlive \
    libgdal-dev \
    ghostscript \
    g++

Codes for Mac OS building enviornment (optional)

The following code is for a local GUI purpose. Hitmap now has been built on the shiny server system. You can skip this step in the later version. You may need to update the Xcode. Go to your Mac OS terminal and input:

xcode-select --install

You’ll then receive: xcode-select: note: install requested for command line developer tools You will be prompted at this point in a window to update Xcode Command Line tools.

You may also need to install the X11.app and tcl/tk support for Mac system:

Example data and source code

The HitMaP comes with a series of maldi-imaging datasets acquired by FT-ICR mass spectromety. With the following code, you can download these raw data set into a local folder.

You can download the example data manually through this link: “https://github.com/MASHUOA/HiTMaP/releases/download/1.0.1/Data.tar.gz

Or download the files in a R console:

if(!require(piggyback)) install.packages("piggyback")
library(piggyback)

#made sure that this folder has enough space
wd="~/expdata/"
dir.create(wd)
setwd(wd)

pb_download("HiTMaP-master.zip", repo = "MASHUOA/HiTMaP", dest = ".",show_progress = F, tag="1.0.1")

pb_download("Data.tar.gz", repo = "MASHUOA/HiTMaP", dest = ".")

untar('Data.tar.gz',exdir =".",  tar="tar")

#unlink('Data.tar.gz')
list.dirs()

The example data contains three folders for three individual IMS datasets, which each contain a configuration file, and the fasta database, respectively: “./Bovinlens_Trypsin_FT” “./MouseBrain_Trypsin_FT” “./Peptide_calibrants_FT”

An Tiny version of data set is also available by using the code below:

if(!require(piggyback)) install.packages("piggyback")
library(piggyback)

#made sure that this folder has enough space
wd="~/expdata/"
dir.create(wd)
setwd(wd)
pb_download("Data_tiny.tar.gz", repo = "MASHUOA/HiTMaP", dest = ".")
untar('Data_tiny.tar.gz',exdir =".",  tar="tar")

#unlink('Data.tar.gz')
list.dirs()

The tiny version dataset was generated from the Bovinlens and MouseBrain original data:

  1. m/z range: 700 - 1400

  2. pixel range:

x \<= 20%, y >= 80% (Bovinlens)

x \<= 30%, y \<= 20% (MouseBrain)

Proteomics identification on maldi-imaging dataset

To perform false-discovery rate controlled peptide and protein annotation, run the following script below in your R session:

#create candidate list
library(HiTMaP)
#set project folder that contains imzML, .ibd and fasta files
#wd=paste0(file.path(path.package(package="HiTMaP")),"/data/")
#set a series of imzML files to be processed
datafile=c("Bovinlens_Trypsin_FT/Bovin_lens.imzML")
wd="~/expdata/"

preprocess = list(force_preprocess=TRUE,
                  use_preprocessRDS=FALSE,
                  smoothSignal=list(method = c("Disable", "gaussian", "sgolay", "ma")[1]),
                  reduceBaseline=list(method = c("Disable", "locmin", "median")[1]),
                  peakPick=list(method=c("diff", "sd", "mad", "quantile", "filter", "cwt")[3]),
                  peakAlign=list(tolerance=5, units="ppm", level=c("local","global")[1], method=c("Enable","Disable")[1]),
                  normalize=list(method=c("Disable","rms","tic","reference")[1], mz=NULL)
                  )

imaging_identification(
#==============Choose the imzml raw data file(s) to process  make sure the fasta file in the same folder
               datafile=paste0(wd,datafile),
               threshold=0.005, 
               ppm=5,
               FDR_cutoff = 0.05,
#==============specify the digestion enzyme specificity
               Digestion_site="trypsin",
#==============specify the range of missed Cleavages
               missedCleavages=0:1,
#==============Set the target fasta file
               Fastadatabase="uniprot-bovin.fasta",
#==============Set the possible adducts and fixed modifications
               adducts=c("M+H"),
               Modifications=list(fixed=NULL,fixmod_position=NULL,variable=NULL,varmod_position=NULL),
#==============The decoy mode: could be one of the "adducts", "elements" or "isotope"
               Decoy_mode = "isotope",
               use_previous_candidates=F,
               output_candidatelist=T,
#==============The pre-processing param
               preprocess=preprocess,
#==============Set the parameters for image segmentation
               spectra_segments_per_file=4,
               Segmentation="spatialKMeans",
               Smooth_range=1,
               Virtual_segmentation=FALSE,
               Virtual_segmentation_rankfile=NULL,
#==============Set the Score method for hi-resolution isotopic pattern matching
               score_method="SQRTP",
               peptide_ID_filter=2,
#==============Summarise the protein and peptide features across the project the result can be found at the summary folder
               Protein_feature_summary=TRUE,
               Peptide_feature_summary=TRUE,
               Region_feature_summary=TRUE,
#==============The parameters for Cluster imaging. Specify the annotations of interest, the program will perform a case-insensitive search on the result file, extract the protein(s) of interest and plot them in the cluster imaging mode
               plot_cluster_image_grid=FALSE,
               ClusterID_colname="Protein",
               componentID_colname="Peptide",
               Protein_desc_of_interest=c("Crystallin","Actin"),
               Rotate_IMG=NULL,
               )

Project folder and result structure

In the above function, you have performed proteomics analysis on the sample data file. It is a tryptic Bovin lens MALDI-imaging file which is acquired on an FT-ICR MS. The function will take the selected data files’ root directory as the project folder. In this example, the project folder will be:

library(HiTMaP)
wd=paste0("D:\\GITHUB LFS\\HiTMaP-Data\\inst","/data/Bovinlens_Trypsin_FT/")
datafile=c("Bovin_lens")

After the whole identification process, you will get two sub-folders within the project folder:

list.dirs(wd, recursive=FALSE)
## [1] "D:\\GITHUB LFS\\HiTMaP-Data\\inst/data/Bovinlens_Trypsin_FT//Bovin_lens ID" 
## [2] "D:\\GITHUB LFS\\HiTMaP-Data\\inst/data/Bovinlens_Trypsin_FT//Summary folder"
  1. The one which has an identical name to an input data file contains the identification result of that input:

    • the protein and peptides list of each segmentation region
    • the PMF matching plot of each segmentation
    • the image that indicates segmentations’ boundary (applies to either K-mean segmentation (powered by Cardinal) or manually defined segmentation)
    • folders of each region contains the detailed identification process, FDR plots and isotopic pattern matching plots
  2. “Summary folder” contains:

    • the identification summary of protein and peptides across all the data
    • the candidate list of all possible proteins and peptides (if use_previous_candidates is set as TRUE)
    • the Cluster imaging files of the protein of interest
    • the database stats result for resolution-based candidates binning (optional)

Identification result visulasation and interpretation

To plot the MALDI-image peptide and protein images, use the following functions:

To check the segmentation result over the sample, you need to navigate to each data file ID folder and find the “spatialKMeans_image_plot.png” (if you are using the spatial K-means method for segmentation.)

library(magick)
p<-image_read(paste0(wd,datafile," ID/spatialKMeans_image_plot.png"))
print(p)
## # A tibble: 1 x 7
##   format width height colorspace matte filesize density
##   <chr>  <int>  <int> <chr>      <lgl>    <int> <chr>  
## 1 PNG     1024   2640 sRGB       FALSE    30726 72x72

The pixels in image data now has been categorized into four regions according to the initial setting of segmentation (spectra_segments_per_file=5). The rainbow shaped bovine lens segmentation image (on the left panel) shows a unique statistical classification based on the mz features of each region (on the right panel).

The mouse brain example segmentation result (spatialKmeans n=9) shown as below:

For further investigation of the segmentation process, you may find a PCA images set in the “Datafile ID” folder. THe PCA images are good summary of features and potential region of interests within a data file. The combination of these PCs of interest will guide you to the insightful tissue structure profile.

The identification will take place on the mean spectra of each region. To check the peptide mass fingerprint (PMF) matching quality, you could locate the PMF spectrum matching plot of each individual region.

library(magick)
p_pmf<-image_read(paste0(wd,datafile," ID/Bovin_lens 3PMF spectrum match.png"))
print(p_pmf)
## # A tibble: 1 x 7
##   format width height colorspace matte filesize density
##   <chr>  <int>  <int> <chr>      <lgl>    <int> <chr>  
## 1 PNG     1980   1080 sRGB       FALSE    17664 72x72

A list of the peptides and proteins annotated within each region has also been created for manual exploration of the results.

peptide_pmf_result<-read.csv(paste0(wd,datafile," ID/Peptide_segment_PMF_RESULT_3.csv"))
head(peptide_pmf_result)
##   Protein       mz Protein_coverage isdecoy       Peptide Modification    pepmz
## 1      48 1300.664       0.06875544       0   HLEQFATEGLR           NA 1299.657
## 2      48 1300.661       0.06875544       0   QYFLDLALSCK           NA 1299.653
## 3      48 1324.643       0.06875544       0   GSKCILYCFYK           NA 1323.636
## 4      53 1328.747       0.05542725       0   FKNINPFPVPR           NA 1327.740
## 5      53 1449.712       0.05542725       0  AVQNFTEYNVHK           NA 1448.705
## 6      53 1605.813       0.05542725       0 AVQNFTEYNVHKR           NA 1604.806
##          formula adduct charge start end pro_end mz_align      Score Rank
## 1   C57H90N17O18    M+H      1   580 590    1149 1300.666  2.4633527    4
## 2 C60H94N13O17S1    M+H      1   744 754    1149 1300.666  2.0216690   10
## 3 C62H94N13O15S2    M+H      1   840 850    1149 1324.647 -0.2644896   32
## 4   C64H98N17O14    M+H      1   207 217     433 1328.747  1.0865820    7
## 5   C65H97N18O20    M+H      1    92 103     433 1449.714  0.7060553   10
## 6  C71H109N22O21    M+H      1    92 104     433 1605.806  2.7178547   11
##   moleculeNames Region Delta_ppm Intensity peptide_count
## 1   HLEQFATEGLR      3 0.9026772 4672324.6             3
## 2   QYFLDLALSCK      3 1.4117311 4672324.6             3
## 3   GSKCILYCFYK      3 1.5164261  145191.4             3
## 4   FKNINPFPVPR      3 0.9094769  191636.4             3
## 5  AVQNFTEYNVHK      3 2.8830137 1275214.1             3
## 6 AVQNFTEYNVHKR      3 1.6464326  558610.4             3
##                                                                                                                                          desc.x
## 1                                  sp|Q29449|AT8A1_BOVIN Probable phospholipid-transporting ATPase IA OS=Bos taurus OX=9913 GN=ATP8A1 PE=1 SV=2
## 2                                  sp|Q29449|AT8A1_BOVIN Probable phospholipid-transporting ATPase IA OS=Bos taurus OX=9913 GN=ATP8A1 PE=1 SV=2
## 3                                  sp|Q29449|AT8A1_BOVIN Probable phospholipid-transporting ATPase IA OS=Bos taurus OX=9913 GN=ATP8A1 PE=1 SV=2
## 4 sp|Q3SX05|ECSIT_BOVIN Evolutionarily conserved signaling intermediate in Toll pathway, mitochondrial OS=Bos taurus OX=9913 GN=ECSIT PE=2 SV=1
## 5 sp|Q3SX05|ECSIT_BOVIN Evolutionarily conserved signaling intermediate in Toll pathway, mitochondrial OS=Bos taurus OX=9913 GN=ECSIT PE=2 SV=1
## 6 sp|Q3SX05|ECSIT_BOVIN Evolutionarily conserved signaling intermediate in Toll pathway, mitochondrial OS=Bos taurus OX=9913 GN=ECSIT PE=2 SV=1
##                                                                                                                                          desc.y
## 1                                  sp|Q29449|AT8A1_BOVIN Probable phospholipid-transporting ATPase IA OS=Bos taurus OX=9913 GN=ATP8A1 PE=1 SV=2
## 2                                  sp|Q29449|AT8A1_BOVIN Probable phospholipid-transporting ATPase IA OS=Bos taurus OX=9913 GN=ATP8A1 PE=1 SV=2
## 3                                  sp|Q29449|AT8A1_BOVIN Probable phospholipid-transporting ATPase IA OS=Bos taurus OX=9913 GN=ATP8A1 PE=1 SV=2
## 4 sp|Q3SX05|ECSIT_BOVIN Evolutionarily conserved signaling intermediate in Toll pathway, mitochondrial OS=Bos taurus OX=9913 GN=ECSIT PE=2 SV=1
## 5 sp|Q3SX05|ECSIT_BOVIN Evolutionarily conserved signaling intermediate in Toll pathway, mitochondrial OS=Bos taurus OX=9913 GN=ECSIT PE=2 SV=1
## 6 sp|Q3SX05|ECSIT_BOVIN Evolutionarily conserved signaling intermediate in Toll pathway, mitochondrial OS=Bos taurus OX=9913 GN=ECSIT PE=2 SV=1
protein_pmf_result<-read.csv(paste0(wd,datafile," ID/Protein_segment_PMF_RESULT_3.csv"))
head(protein_pmf_result)
##   Protein   Proscore isdecoy Intensity     Score peptide_count Protein_coverage
## 1   10134 0.13943597       0 2873903.1 1.9269417             3       0.06715328
## 2   10204 0.13654123       0  380571.3 0.7940642             3       0.18468468
## 3   10370 0.20365140       0 1877250.1 2.0776861             4       0.09364548
## 4   10659 0.11239668       0  327352.4 0.7448240             3       0.16400000
## 5   10888 0.07975644       0  532832.0 1.2420183             3       0.06720978
## 6   11270 0.10687770       0 2944154.2 1.3292158             3       0.07449857
##   Intensity_norm
## 1      1.0775539
## 2      0.9310593
## 3      1.0466962
## 4      0.9201443
## 5      0.9554442
## 6      1.0793038
##                                                                                                                          desc
## 1                   tr|G3N2M8|G3N2M8_BOVIN Sterile alpha motif domain containing 15 OS=Bos taurus OX=9913 GN=SAMD15 PE=4 SV=2
## 2                                      tr|A0A3Q1LYB6|A0A3Q1LYB6_BOVIN Uncharacterized protein OS=Bos taurus OX=9913 PE=4 SV=1
## 3             tr|E1B9U7|E1B9U7_BOVIN Polypeptide N-acetylgalactosaminyltransferase OS=Bos taurus OX=9913 GN=GALNT17 PE=3 SV=3
## 4 tr|A0A3Q1M1B1|A0A3Q1M1B1_BOVIN Phosphatidylinositol transfer protein beta isoform OS=Bos taurus OX=9913 GN=PITPNB PE=4 SV=1
## 5                                  tr|F1MMD4|F1MMD4_BOVIN Matrix metallopeptidase 11 OS=Bos taurus OX=9913 GN=MMP11 PE=3 SV=2
## 6                         tr|F6RR01|F6RR01_BOVIN Ribosome production factor 1 homolog OS=Bos taurus OX=9913 GN=RPF1 PE=4 SV=1

Scoring system for protein and peptide

Score in peptide result table shows the isotopic pattern matching score of the peptide (Pepscore). In Protein result table, it shows the protein score (Proscore). The ‘Pepscore’ consist of two parts: Intensity_Score and Mass_error_Score:

Proscore in the protein result table shows the overall estimation of the protein identification Accuracy.

A Peptide_region_file.csv has also been created to summarise all the IDs in this data file:

Identification_summary_table<-read.csv(paste0(wd,datafile," ID/Peptide_region_file.csv"))
head(Identification_summary_table)
##   Protein        mz Protein_coverage isdecoy              Peptide Modification
## 1      24 1143.5793       0.06119704       0         GFPGQDGLAGPK           NA
## 2      24 1684.8878       0.06119704       0   DGANGIPGPIGPPGPRGR           NA
## 3      24  742.3478       0.06119704       0             GDSGPPGR           NA
## 4      24 1693.8214       0.06119704       0      LLSTEGSQNITYHCK           NA
## 5      24 1881.9276       0.06119704       0 GQPGVMGFPGPKGANGEPGK           NA
## 6      48 1216.7008       0.03481288       0          ASTSVQNRLLK           NA
##       pepmz         formula adduct charge start  end pro_end  mz_align
## 1 1142.5720    C51H79N14O16    M+H      1   516  527    1487 1143.5828
## 2 1683.8805   C72H118N25O22    M+H      1  1175 1192    1487 1684.8830
## 3  741.3406    C29H48N11O12    M+H      1   933  940    1487  742.3504
## 4 1692.8141 C72H117N20O25S1    M+H      1  1380 1394    1487 1693.8197
## 5 1880.9203 C82H129N24O25S1    M+H      1   597  616    1487 1881.9268
## 6 1215.6935    C51H94N17O17    M+H      1   614  624    1149 1216.7047
##       Score Rank        moleculeNames Region Delta_ppm Intensity peptide_count
## 1 1.4443497    2         GFPGQDGLAGPK      2 1.3471596  250698.3             5
## 2 1.9337304    2   DGANGIPGPIGPPGPRGR      2 1.5937657 2696717.3             5
## 3 1.2698949    1             GDSGPPGR      2 0.1407633  190469.7             5
## 4 1.3660521    3      LLSTEGSQNITYHCK      2 2.2329023  368927.9             5
## 5 0.5868561   17 GQPGVMGFPGPKGANGEPGK      2 3.0817671  974427.3             5
## 6 1.9039495    1          ASTSVQNRLLK      2 1.8837090 2036000.7             1
##                                                                                                         desc.x
## 1                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 2                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 3                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 4                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 5                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 6 sp|Q29449|AT8A1_BOVIN Probable phospholipid-transporting ATPase IA OS=Bos taurus OX=9913 GN=ATP8A1 PE=1 SV=2
##                                                                                                         desc.y
## 1                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 2                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 3                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 4                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 5                   sp|P02459|CO2A1_BOVIN Collagen alpha-1(II) chain OS=Bos taurus OX=9913 GN=COL2A1 PE=1 SV=4
## 6 sp|Q29449|AT8A1_BOVIN Probable phospholipid-transporting ATPase IA OS=Bos taurus OX=9913 GN=ATP8A1 PE=1 SV=2

The details of protein/peptide identification process has been save to the folder named by the segmentation:

list.dirs(paste0(wd,datafile," ID/"), recursive=FALSE)
## [1] "D:\\GITHUB LFS\\HiTMaP-Data\\inst/data/Bovinlens_Trypsin_FT/Bovin_lens ID//1"
## [2] "D:\\GITHUB LFS\\HiTMaP-Data\\inst/data/Bovinlens_Trypsin_FT/Bovin_lens ID//2"
## [3] "D:\\GITHUB LFS\\HiTMaP-Data\\inst/data/Bovinlens_Trypsin_FT/Bovin_lens ID//3"
## [4] "D:\\GITHUB LFS\\HiTMaP-Data\\inst/data/Bovinlens_Trypsin_FT/Bovin_lens ID//4"

In the identification details folder, you will find a series of FDR file and plots to demonstrate the FDR model and score cutoff threshold:

dir(paste0(wd,datafile," ID/1/"), recursive=FALSE)
##  [1] "FDR.CSV"                                        
##  [2] "FDR.png"                                        
##  [3] "Matching_Score_vs_mz_target-decoy.png"          
##  [4] "Peptide_1st_ID.csv"                             
##  [5] "Peptide_1st_ID_score_rank_SQRTP.csv"            
##  [6] "Peptide_2nd_ID_score_rankSQRTP_Rank_above_3.csv"
##  [7] "Peptide_Score_histogram_target-decoy.png"       
##  [8] "ppm"                                            
##  [9] "PROTEIN_FDR.CSV"                                
## [10] "Protein_FDR.png"                                
## [11] "Protein_ID_score_rank_SQRTP.csv"                
## [12] "PROTEIN_Score_histogram.png"                    
## [13] "Spectrum.csv"                                   
## [14] "unique_peptide_ranking_vs_mz_feature.png"

In this folder, you will find the FDR plots for protein and peptide annotation. The software will take the proscore and its FDR model to trim the final identification results. The unique_peptide_ranking_vs_mz_feature.png is a plot that could tell you the number of peptide candidates that have been matched to the mz features in the first round run. You can also access the peptide spectrum match (first MS dimension) data via the “/ppm” subfolder.

library(magick)
p_FDR_peptide<-image_read(paste0(wd,datafile," ID/3/FDR.png"))
p_FDR_protein<-image_read(paste0(wd,datafile," ID/3/protein_FDR.png"))
p_FDR_peptide_his<-image_read(paste0(wd,datafile," ID/3/Peptide_Score_histogram_target-decoy.png"))
p_FDR_protein_his<-image_read(paste0(wd,datafile," ID/3/PROTEIN_Score_histogram.png"))
p_combined<-image_append(c(p_FDR_peptide,p_FDR_peptide_his,p_FDR_protein,p_FDR_protein_his))
print(p_combined)
## # A tibble: 1 x 7
##   format width height colorspace matte filesize density
##   <chr>  <int>  <int> <chr>      <lgl>    <int> <chr>  
## 1 PNG     1920    480 sRGB       FALSE        0 72x72

You will also find a Matching_Score_vs_mz plot for further investigation on peptide matching quality.

library(magick)
#plot Matching_Score_vs_mz
p_Matching_Score_vs_mz<-image_read(paste0(wd,datafile," ID/3/Matching_Score_vs_mz_target-decoy.png"))
print(p_Matching_Score_vs_mz)
## # A tibble: 1 x 7
##   format width height colorspace matte filesize density
##   <chr>  <int>  <int> <chr>      <lgl>    <int> <chr>  
## 1 PNG      480    480 sRGB       FALSE    47438 72x72

Identification summary and cluster imaging

In the project summary folder, you will find four files and a sub-folder:

wd_sum=paste(wd,"/Summary folder", sep="")
dir(wd_sum)
## [1] "candidatelist.csv"                "cluster Ion images"              
## [3] "Peptide_Summary.csv"              "Protein_feature_list_trimmed.csv"
## [5] "protein_index.csv"                "Protein_Summary.csv"

“candidatelist.csv” and “protein_index.csv” contains the candidates used for this analysis. They are output after the candidate processing while output_candidatelist set as TRUE, and can be used repeatedly while use_previous_candidates set as TRUE.

We have implemented a functionality to perform additional statistical analysis around the number of enzymatically generated peptides derived from a given proteome database. If the user sets the argument ‘Database_stats’ to TRUE in the main workflow, the function will be called. Briefly, the function will list all of the m/z’s of a unique formulae from the peptide candidate pool within a given m/z range. The m/z’s will then be binned using three tolerance window: 1 ppm, 2 ppm and 5 ppm. A plot showing the number of unique formulae vs. m/z bins will be generated and exported to the summary folder (DB_stats_mz_bin).

Proteome database stats

“Peptide_Summary.csv” and “Protein_Summary.csv” contains the table of the project identification summary. You could set the plot_cluster_image_grid as TRUE to enable the cluster imaging function. Please be noted that you could indicate Rotate_IMG with a CSV file path that indicates the rotation degree of image files.

Note: 90°, 180° and 270° are recommended for image rotation. You may find an example CSV file in the expdata/MouseBrain_Trypsin_FT/file_rotationbk.csv.

library(dplyr)
Protein_desc_of_interest<-c("Crystallin","Actin")
Protein_Summary_tb<-read.csv(paste(wd,"/Summary folder","/Protein_Summary.csv", sep=""),stringsAsFactors = F)

Finally, you are able visualize the annotated proteins and their associated peptide distributions via the cluster image rendering function.

Bovin lens cluster image rendering

vimentin:

β-crystallin:

α-crystallin:

Secernin 1

CX6A1 cytochrome coxidase subunit 6A1

Myelin basic protein

Pixel level Proteomics data export

Details of parameter setting

Modification

You can choose one or a list of modifications from the unimod modification list. Peptide_modification function is used to load/rebuild the modification database into the global enviornment of R. It will be called automatically in the identification work flow. you can use the code_name or record_id to refer the modification (see example data “peptide calibrants” to find more details). The pipeline will select the non-hidden modifications.

HiTMaP:::Peptide_modification(retrive_ID=NULL,update_unimod=F)
modification_list<-merge(unimod.df$modifications,unimod.df$specificity,by.x=c("record_id"),by.y=c("mod_key"),all.x=T)
head(modification_list['&'(modification_list$code_name=="Phospho",modification_list$hidden!=1),c("code_name","record_id","composition","mono_mass","position_key","one_letter")])
##      code_name record_id composition mono_mass position_key one_letter
## 1615   Phospho        21    H O(3) P 79.966331            2          Y
## 1618   Phospho        21    H O(3) P 79.966331            2          T
## 1620   Phospho        21    H O(3) P 79.966331            2          S
head(modification_list['&'(modification_list$code_name=="Amide",modification_list$hidden!=1),c("code_name","record_id","composition","mono_mass","position_key","one_letter")])
##      code_name record_id composition mono_mass position_key one_letter
## 1552     Amide         2   H N O(-1) -0.984016            6     C-term
## 1553     Amide         2   H N O(-1) -0.984016            4     C-term
head(modification_list['&'(stringr::str_detect(modification_list$code_name,"Ca"),modification_list$hidden!=1),c("code_name","record_id","composition","mono_mass","position_key","one_letter")])
##            code_name record_id    composition mono_mass position_key one_letter
## 1946 Carbamidomethyl         4  H(3) C(2) N O 57.021464            2          C
## 1949 Carbamidomethyl         4  H(3) C(2) N O 57.021464            3     N-term
## 2119        Carbamyl         5        H C N O 43.005814            3     N-term
## 2121        Carbamyl         5        H C N O 43.005814            2          K
## 2271   Carboxymethyl         6 H(2) C(2) O(2) 58.005479            2          C

If a modification occurs on a particular site, you will also need to specify the position of a modifications.

unimod.df[["positions"]]
##         position record_id
## 1              -         1
## 2       Anywhere         2
## 3     Any N-term         3
## 4     Any C-term         4
## 5 Protein N-term         5
## 6 Protein C-term         6

Amino acid substitution

You can set the Substitute_AA to make the uncommon amino acid available to the workflow: Substitute_AA=list(AA=c(“X”),AA_new_formula=c(“C5H5NO2”),Formula_with_water=c(FALSE))

Digestion site and enzyme

The Digestion_site allows you to specify a list of pre-defined enzyme and customized digestion rules in regular expression format. You can either use the enzyme name, customized cleavage rule or combination of them to get the enzymatics peptides list.

Cleavage_rules<-Cleavage_rules_fun()
Cleavage_df<-data.frame(Enzyme=names(Cleavage_rules),Cleavage_rules=unname(Cleavage_rules),stringsAsFactors = F)
library(gridExtra)
grid.ftable(Cleavage_df, gp = gpar(fontsize=9,fill = rep(c("grey90", "grey95"))))

Imaging-MS data preprocessing

preprocess\$mz_bin_list is an argument for costumized peak-picking and mz bining purpose. If it is not NULL, the workflow will bypass signal smooth, noise reduction, and peakpicking steps. User need to give a numeric vector as the mz input to this argument. The workflow will first filter the vector with the given ppm tolerance to ensure there’s no overlapped mz bins (mz +/- ppm tolerance). Then, a m/z binning procedure will be conducted to the image data to produce a peak-picked dataset (the peak bin width will be the ppm tolerance).

If user uses a processed IMS data that contains the centroid feature value (e.g. exported from scils lab with feature list reduced data). User will still be safe to use this mz_bin_list in order to mount the centroid data properly. In this case, the ppm tolerance will only applied to the following annotation procedure.

normalize=list(method=c(“Disable”,“rms”,“tic”,“reference”)[1],mz=1) the current IMS normalization is done on pixel-to-pixel level, which will affect the feature distribution in some tissue. We use “Disable” in the example dataset to minimize the required RAM space and working time. The step may result in a big RAM usage on some IMS data. If the error message mentioned a “vector allocation” issue, Please consider to disable the normalization.

Example workflow command

Below is a list of commands including the parameters for the example data sets.

Peptide calibrant

#peptide calibrant
library(HiTMaP)
datafile=c("Peptide_calibrants_FT/trypsin_non-decell_w.calibrant_FTICR")
wd="~/expdata/"

# Calibrants dataset analysis with modification
imaging_identification(datafile=paste0(wd,datafile),
  Digestion_site="trypsin",
  Fastadatabase="uniprot_cali.fasta",
  output_candidatelist=T,
  plot_matching_score=T,
  spectra_segments_per_file=1,
  use_previous_candidates=F,
  peptide_ID_filter=1,ppm=5,missedCleavages=0:5,
  Modifications=list(fixed=NULL,fixmod_position=NULL,variable=c("Amide"),varmod_position=c(6)),
  FDR_cutoff=0.1,
  Substitute_AA=list(AA=c("X"),AA_new_formula=c("C5H5NO2"),Formula_with_water=c(FALSE)))

# Calibrants dataset analysis with no modification
imaging_identification(datafile=paste0(wd,datafile),
  Digestion_site="trypsin",
  Fastadatabase="uniprot_cali.fasta",
  output_candidatelist=T,
  plot_matching_score=T,
  spectra_segments_per_file=1,
  use_previous_candidates=T,
  peptide_ID_filter=1,ppm=5,missedCleavages=0:5,
  FDR_cutoff=0.1)

library(HiTMaP)
datafile=c("Peptide_calibrants_FT/trypsin_non-decell_w.calibrant_FTICR")
wd="~/expdata/"
# Calibrants dataset analysis with modification 
imaging_identification(datafile=paste0(wd,datafile),
  Digestion_site="trypsin",
  Fastadatabase="calibrants.fasta",
  output_candidatelist=T,
  plot_matching_score=T,
  spectra_segments_per_file=1,
  use_previous_candidates=F,
  peptide_ID_filter=1,ppm=5,missedCleavages=0:5,
  Modifications=list(fixed=NULL,fixmod_position=NULL,variable=c("Amide"),varmod_position=c(6)),
  FDR_cutoff=0.1,
  Substitute_AA=list(AA=c("X"),AA_new_formula=c("C5H5NO2"),Formula_with_water=c(FALSE)),Thread = 1)

Bovine lens

library(HiTMaP)
datafile=c("Bovinlens_Trypsin_FT/Bovin_lens.imzML")
wd="~/expdata/"

# Data pre-processing and proteomics annotation
library(HiTMaP)
imaging_identification(datafile=paste0(wd,datafile),Digestion_site="trypsin",
                       Fastadatabase="uniprot-bovin.fasta",output_candidatelist=T,
                       preprocess=list(force_preprocess=TRUE,
                               use_preprocessRDS=TRUE,
                               smoothSignal=list(method="Disable"),
                               reduceBaseline=list(method="Disable"),
                               peakPick=list(method="adaptive"),
                               peakAlign=list(tolerance=5, units="ppm"),
                               normalize=list(method=c("Disable","rms","tic","reference")[1],mz=1)),
                       spectra_segments_per_file=4,use_previous_candidates=F,ppm=5,FDR_cutoff = 0.05,IMS_analysis=T,
                       Rotate_IMG="file_rotationbk.csv",plot_cluster_image_grid=F)

datafile=c("Bovinlens_Trypsin_FT/Bovin_lens.imzML")
wd="~/expdata/"
library(HiTMaP)
imaging_identification(datafile=paste0(wd,datafile),Digestion_site="trypsin",
                       Fastadatabase="uniprot-bovin.fasta",output_candidatelist=T,use_previous_candidates=T,
                       preprocess=list(force_preprocess=F,
                               use_preprocessRDS=TRUE,
                               smoothSignal=list(method="Disable"),
                               reduceBaseline=list(method="Disable"),
                               peakPick=list(method="Default"),
                               peakAlign=list(tolerance=5, units="ppm"),
                               normalize=list(method=c("Disable","rms","tic","reference")[1],mz=1)),
                       spectra_segments_per_file=4,ppm=5,FDR_cutoff = 0.05,IMS_analysis=T,
                       Rotate_IMG="file_rotationbk.csv",plot_cluster_image_grid=F)

# Re-analysis and cluster image rendering

library(HiTMaP)
datafile=c("Bovinlens_Trypsin_FT/Bovin_lens.imzML")
wd="~/expdata/"
imaging_identification(datafile=paste0(wd,datafile),Digestion_site="trypsin",
                       Fastadatabase="uniprot-bovin.fasta",
                       use_previous_candidates=T,ppm=5,IMS_analysis=F,
                       plot_cluster_image_grid=T,
                       export_Header_table=T, 
                       img_brightness=250, 
                       plot_cluster_image_overwrite=T,
                       cluster_rds_path = "/Bovin_lens ID/preprocessed_imdata.RDS",pixel_size_um = 150,
                       Plot_score_abs_cutoff=-0.1,
                       remove_score_outlier=T,
                       Protein_desc_of_interest=c("Crystallin","Phakinin","Filensin","Actin","Vimentin","Cortactin","Visinin","Arpin","Tropomyosin","Myosin Light Chain 3","Kinesin Family Member 14","Dynein Regulatory Complex","Ankyrin Repeat Domain 45"))

# Re-analysis and cluster image rendering using color scale

library(HiTMaP)
datafile=c("Bovinlens_Trypsin_FT/Bovin_lens.imzML")
wd="~/expdata/"
imaging_identification(datafile=paste0(wd,datafile),Digestion_site="trypsin",
                       Fastadatabase="uniprot-bovin.fasta",
                       use_previous_candidates=T,ppm=5,IMS_analysis=F,
                       plot_cluster_image_grid=T,
                       export_Header_table=T, 
                       img_brightness=250, 
                       plot_cluster_image_overwrite=T,
                       cluster_rds_path = "/Bovin_lens ID/preprocessed_imdata.RDS",pixel_size_um = 150,
                       Plot_score_abs_cutoff=-0.1,
                       remove_score_outlier=T,cluster_color_scale="fleximaging",
                       Protein_desc_of_interest=c("Crystallin","Phakinin","Filensin","Actin","Vimentin","Cortactin","Visinin","Arpin","Tropomyosin","Myosin Light Chain 3","Kinesin Family Member 14","Dynein Regulatory Complex","Ankyrin Repeat Domain 45"))

Mouse brain

library(HiTMaP)
datafile=c("MouseBrain_Trypsin_FT/Mouse_brain.imzML")
wd="~/expdata/"
preprocess = list(force_preprocess=TRUE,
                  use_preprocessRDS=FALSE,
                  smoothSignal=list(method = c("Disable", "gaussian", "sgolay", "ma")[1]),
                  reduceBaseline=list(method = c("Disable", "locmin", "median")[1]),
                  peakPick=list(method=c("diff", "sd", "mad", "quantile", "filter", "cwt")[3]),
                  peakAlign=list(tolerance=5, units="ppm", level=c("local","global")[1], method=c("Enable","Disable")[1]),
                  normalize=list(method=c("Disable","rms","tic","reference")[1], mz=NULL)
                  )

# Data pre-processing and proteomics annotation
library(HiTMaP)
imaging_identification(datafile=paste0(wd,datafile),Digestion_site="trypsin",
                       Fastadatabase="uniprot_mouse_20210107.fasta",output_candidatelist=T,
                       preprocess=preprocess,
                       spectra_segments_per_file=9,use_previous_candidates=F,ppm=10,FDR_cutoff = 0.05,IMS_analysis=T,
                       Rotate_IMG="file_rotationbk.csv",
                       mzrange = c(500,4000),plot_cluster_image_grid=F)

# Re-analysis and cluster image rendering
library(HiTMaP)
datafile=c("MouseBrain_Trypsin_FT/Mouse_brain.imzML")
wd="~/expdata/"
imaging_identification(datafile=paste0(wd,datafile),Digestion_site="trypsin",
                       Fastadatabase="uniprot_mouse_20210107.fasta",
                       preprocess=list(force_preprocess=FALSE),
                       spectra_segments_per_file=9,use_previous_candidates=T,ppm=10,FDR_cutoff = 0.05,IMS_analysis=F,
                       mzrange = c(500,4000),plot_cluster_image_grid=T,
                       img_brightness=250, plot_cluster_image_overwrite=T,
                       cluster_rds_path = "/Mouse_brain ID/preprocessed_imdata.RDS",
                       pixel_size_um = 50,
                       Plot_score_abs_cutoff=-0.1,
                       remove_score_outlier=T,
                       Protein_desc_of_interest=c("Secernin","GN=MBP","Cytochrome"))

library(HiTMaP)
datafile=c("MouseBrain_Trypsin_FT_200brit_man_seg/Mouse_brain.imzML")
wd="~/expdata/"
imaging_identification(datafile=paste0(wd,datafile),Digestion_site="trypsin",
                       Fastadatabase="uniprot_mouse_20210107.fasta",
                       preprocess=list(force_preprocess=FALSE),
                       spectra_segments_per_file=9,use_previous_candidates=T,ppm=10,FDR_cutoff = 0.05,IMS_analysis=F,
                       mzrange = c(500,4000),plot_cluster_image_grid=T,
                       img_brightness=250, plot_cluster_image_overwrite=T,
                       cluster_rds_path = "/Mouse_brain ID/preprocessed_imdata.RDS",
                       pixel_size_um = 50,
                       Plot_score_abs_cutoff=-0.1,
                       remove_score_outlier=T,
                       Protein_desc_of_interest=c("GTR9"))

Cite this project

This study has been accepted by Nature Communications: <DOI:10.1038/s41467-021-23461-w> “Automated annotation and visualisation of high-resolution spatial proteomic mass spectrometry imaging data using HIT-MAP” online on the 28th May 2021.

Session information

sessionInfo()
## R version 4.4.2 (2024-10-31 ucrt)
## Platform: x86_64-w64-mingw32/x64
## Running under: Windows 10 x64 (build 19045)
## 
## Matrix products: default
## 
## 
## locale:
## [1] LC_COLLATE=English_Australia.utf8  LC_CTYPE=English_Australia.utf8   
## [3] LC_MONETARY=English_Australia.utf8 LC_NUMERIC=C                      
## [5] LC_TIME=English_Australia.utf8    
## 
## time zone: Pacific/Auckland
## tzcode source: internal
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] gridExtra_2.3 XML_3.99-0.17 protViz_0.7.9 dplyr_1.1.4   magick_2.8.5 
## [6] HiTMaP_1.0.0 
## 
## loaded via a namespace (and not attached):
##  [1] vctrs_0.6.5       cli_3.6.3         knitr_1.49        rlang_1.1.4      
##  [5] xfun_0.49         generics_0.1.3    glue_1.8.0        htmltools_0.5.8.1
##  [9] fansi_1.0.6       rmarkdown_2.29    grid_4.4.2        evaluate_1.0.1   
## [13] tibble_3.2.1      fastmap_1.2.0     yaml_2.3.10       lifecycle_1.0.4  
## [17] compiler_4.4.2    codetools_0.2-20  Rcpp_1.0.13-1     pkgconfig_2.0.3  
## [21] rstudioapi_0.17.1 digest_0.6.37     R6_2.5.1          tidyselect_1.2.1 
## [25] utf8_1.2.4        pillar_1.9.0      magrittr_2.0.3    gtable_0.3.6     
## [29] tools_4.4.2

End of the tutorial, Enjoy~

References

R Packages used in this project: