TL;DR

suppressMessages(require(netDx))
suppressMessages(require(GenomicRanges))
suppressMessages(require(biomaRt)) # for fetching gene coordinates

# read patient CNVs
phenoFile <- sprintf("%s/extdata/AGP1_CNV.txt",path.package("netDx"))
pheno   <- read.delim(phenoFile,sep="\t",header=TRUE,as.is=TRUE)
# sample metadata table must have ID and STATUS columns
colnames(pheno)[1] <- "ID"

# create GRanges object. 
# Must have ID and LOCUS_NAMES in metadata
cnv_GR    <- GRanges(pheno$seqnames,
        IRanges(pheno$start,pheno$end),
        ID=pheno$ID,LOCUS_NAMES=pheno$Gene_symbols)
pheno <- pheno[!duplicated(pheno$ID),]

pathFile <- fetchPathwayDefinitions(
    "February",2018,verbose=TRUE)
pathwayList <- readPathways(pathFile)

# get gene coordinates, use hg18
ensembl <- suppressMessages(suppressWarnings(
   useMart("ENSEMBL_MART_ENSEMBL",
    dataset="hsapiens_gene_ensembl",
    host="may2009.archive.ensembl.org",
    path="/biomart/martservice",archive=FALSE)
))
genes <- getBM(attributes=c("chromosome_name",
        "start_position",
        "end_position",
        "hgnc_symbol"),
    mart=ensembl)
genes <- genes[which(genes[,4]!=""),]
gene_GR     <- GRanges(genes[,1],
        IRanges(genes[,2],genes[,3]),
    name=genes[,4]
)

# create GRangesList of pathway ranges
path_GRList <- mapNamedRangesToSets(gene_GR,pathwayList)

outDir <- sprintf("%s/200129_threeWay",tempdir())
if (file.exists(outDir)) unlink(outDir,recursive=TRUE); dir.create(outDir)

predictClass    <- "case"
out <- suppressMessages(
   buildPredictor_sparseGenetic(
            pheno, cnv_GR, predictClass,
      path_GRList,outDir,
      numSplits=3L, featScoreMax=3L,
      enrichLabels=TRUE,numPermsEnrich=20L,
      numCores=2L)
)

# plot ROC curve. Note that the denominator only includes
# patients with events in networks that are label-enriched
dat <- out$performance_denEnrichedNets
plot(0,0,type="n",xlim=c(0,100),ylim=c(0,100),
    las=1, xlab="False Positive Rate (%)", 
    ylab="True Positive Rate (%)",
    bty='n',cex.axis=1.5,cex.lab=1.3,
    main="ROC curve - Patients in label-enriched pathways")
points(dat$other_pct,dat$pred_pct,
      col="red",type="o",pch=16,cex=1.3,lwd=2)

# calculate AUROC and AUPR
tmp <- data.frame(  
    score=dat$score,
    tp=dat$pred_ol,fp=dat$other_ol,
    # tn: "-" that were correctly not called
    tn=dat$other_tot - dat$other_ol,
    # fn: "+" that were not called 
    fn=dat$pred_tot - dat$pred_ol) 

stats <- netDx::perfCalc(tmp)
tmp <- stats$stats
message(sprintf("PRAUC = %1.2f\n", stats$prauc))
message(sprintf("ROCAUC = %1.2f\n", stats$auc))

# examine pathway-level scores; these are 
# cumulative across the splits - here, each of three
# splits has a max feature score of three, so
# a feature can score a max of (3 + 3 + 3) = 9.
print(head(out$cumulativeFeatScores))

Introduction

netDx natively handles missing data, making it suitable to build predictors with sparse genetic data such as somatic DNA mutations, frequently seen in cancer, and from DNA Copy Number Variations (CNV). This example demonstrates how to use netDx to build a predictor from sparse genetic data. Here we build a case/control classifier for Autism Spectrum Disorder (ASD) diagnosis, starting from rare CNVs. The data is from [Pinto et al. (2014) AJHG 94:677)(https://pubmed.ncbi.nlm.nih.gov/24768552-convergence-of-genes-and-cellular-pathways-dysregulated-in-autism-spectrum-disorders/) .

Design and Adapting the Algorithm for Sparse Event Data

In this design, we group CNVs by pathways. The logic behind the grouping is prior evidence showing that genetic events in diseases tend to converge on cellular processes of relevance to the pathophysiology of the disease. For example, see the Pinto et al. paper referenced in the previous section.

Label enrichment

In this design, similarity is defined as a binary function, a strategy that has advantages and drawbacks. In plain terms, if two patients share a mutation in a pathway, their similarity for that pathway is 1.0 ; otherwise it is zero. This binary definition, while conceptually intuitive, increases the false positive rate in the netDx feature selection step. That is, networks with even a single case sample will get a high feature score, regardless of whether that network is enriched for case samples.

To counter this problem, we introduce a label-enrichment step in the feature selection. A bias measure is first computed for each network, such that a network with only cases has +1; one with only controls has a score of -1; and one with an equal number of both has a score of zero. Label-enrichment compares the bias in each real network, to the bias in that network in label-permuted data. It then assigns an empirical p-value for the proportion of times a label-permuted network has a bias as high as the real network. Only networks with a p-value below a user-assigned threshold pass label-enrichment, and feature selection is limited to these networks. In netDx, label-enrichment is enabled by setting enrichLabels=TRUE in the call to buildPredictor_sparseGenetic().

Cumulative feature scoring

The other difference between this design and those with non-sparse data, is the method of feature scoring. The user specifies a parameter which indicates the number of times to split the data and run feature selection. The algorithm then runs feature selection numSplits times, each time leaving 1/numSplits of the samples out. In each split, features are scored between 0 and featScoreMax, using the same approach as is used for continuous-valued input. Feature scores are then added across the splits so that a feature can score as high as numSplits*featScoreMax.

Evaluating model performance

For a given cutoff for features, a patient is called a “case” if they have a genetic event in pathways that pass feature selection at that cutoff; otherwise, at that cutoff, they are labelled a “control”. These calls are used to generate the false positive and true positive rates across the various cutoffs, which ultimately generates an ROC curve.

Setup

suppressMessages(require(netDx))
suppressMessages(require(GenomicRanges))
suppressMessages(require(biomaRt)) # for fetching gene coordinates

Data

CNV coordinates are read in, and converted into a GRanges object. As always, the sample metadata table, here the pheno object, must have ID and STATUS columns.

outDir <- sprintf("%s/200129_threeWay",tempdir())
if (file.exists(outDir)) unlink(outDir,recursive=TRUE); 
dir.create(outDir)

cat("* Setting up sample metadata\n")
## * Setting up sample metadata
phenoFile <- sprintf("%s/extdata/AGP1_CNV.txt",path.package("netDx"))
pheno   <- read.delim(phenoFile,sep="\t",header=TRUE,as.is=TRUE)
colnames(pheno)[1] <- "ID"
head(pheno)
##        ID seqnames    start      end    Gene_symbols Pathogenic STATUS
## 3  1020_4     chr3  4110452  4145874                         no   case
## 4  1030_3    chr10 56265896 56361311                         no   case
## 5  1030_3     chr7 64316996 64593616 ZNF92,LOC441242         no   case
## 7  1045_3     chr3 83206919 83239473                         no   case
## 11 1050_3     chr6 57021412 57062509        KIAA1586         no   case
## 16 1116_4     chr1 30334653 30951250                         no   case
cnv_GR    <- GRanges(pheno$seqnames,IRanges(pheno$start,pheno$end),
                        ID=pheno$ID,LOCUS_NAMES=pheno$Gene_symbols)
pheno <- pheno[!duplicated(pheno$ID),]

Group CNVs by pathways

The fetchPathwayDefinitions() function downloads pathway definitions from baderlab.org but users may provide custom .gmt files as well. In the example below, gene coordinates for the hg18 genome build are also automatically fetched from BioMart, using the biomaRt package, and converted to a GRanges object. The function mapNamedRangesToSets() is used to group this GRanges object into pathway-level sets.

pathFile <- fetchPathwayDefinitions("February",2018,verbose=TRUE)
## Fetching http://download.baderlab.org/EM_Genesets/February_01_2018/Human/symbol/Human_AllPathways_February_01_2018_symbol.gmt
## Downloading example pathway definitions (only required once)
pathwayList <- readPathways(pathFile)
## ---------------------------------------
## File: 35f2679114c4_Human_AllPathways_February_01_2018_symbol.gmt
## Read 3199 pathways in total, internal list has 3163 entries
##  FILTER: sets with num genes in [10, 200]
##    => 1044 pathways excluded
##    => 2119 left
# get gene coordinates, use hg18
ensembl <- suppressMessages(suppressWarnings(
   useMart("ENSEMBL_MART_ENSEMBL",
    dataset="hsapiens_gene_ensembl",
    host="may2009.archive.ensembl.org",
    path="/biomart/martservice",archive=FALSE)
))
genes <- getBM(attributes=c("chromosome_name",
        "start_position",
        "end_position",
        "hgnc_symbol"),
    mart=ensembl)
genes <- genes[which(genes[,4]!=""),]
gene_GR     <- GRanges(genes[,1],IRanges(genes[,2],genes[,3]),
   name=genes[,4])

Group gene extents into pathway-based sets, which effectively creates grouping rules for netDx. The function mapNamedRangesToSets() does this grouping, generating a GRangesList object.

path_GRList <- mapNamedRangesToSets(gene_GR,pathwayList)

Run predictor

Once the phenotype matrix and grouping rules are set up, the predictor is called using buildPredictor_sparseGenetic(). Note that unlike with non-sparse data, the user does not provide a custom similarity function in this application; currently, the only option available is the binary similarity defined above. As discussed above, setting enrichLabels=TRUE to enable label-enrichment is highly recommended to reduce false positive rate.

predictClass    <- "case"
out <- suppressMessages(
   buildPredictor_sparseGenetic(pheno, cnv_GR, predictClass,
                             path_GRList,outDir,
                             numSplits=3L, featScoreMax=3L,
                             enrichLabels=TRUE,numPermsEnrich=20L,
                             numCores=2L)
)
## Time difference of 12.97821 secs
## Time difference of 10.28171 secs
##          TT_STATUS
## STATUS    TEST TRAIN
##   case     188   376
##   control  208   418
## [1] 794
##    user  system elapsed 
##   0.815   0.209  13.267 
##    Min. 1st Qu.  Median    Mean 3rd Qu.    Max. 
## -1.0000 -0.6667  0.2000  0.1196  1.0000  1.0000 
## [1] 440
## Time difference of 8.45909 secs
##          TT_STATUS
## STATUS    TEST TRAIN
##   case     188   376
##   control  208   418
## [1] 794
##    user  system elapsed 
##   0.644   0.160  11.289 
##    Min. 1st Qu.  Median    Mean 3rd Qu.    Max. 
## -1.0000 -0.6667  0.1505  0.1000  1.0000  1.0000 
## [1] 405
## Time difference of 6.626661 secs
##          TT_STATUS
## STATUS    TEST TRAIN
##   case     188   376
##   control  210   416
## [1] 792
##    user  system elapsed 
##   0.792   0.232  11.078 
##    Min. 1st Qu.  Median    Mean 3rd Qu.    Max. 
## -1.0000 -0.6000  0.3333  0.2174  1.0000  1.0000 
## [1] 378
## Time difference of 6.700725 secs

Plot results

Feature selection identifies pathways that are consistently enriched for the label of interest; here, “case” status. From the diagnostic point of view, a patient with a genetic event in a selected feature - here, a CNV in a feature-selected pathway - is labelled a “case”. “True positives” are therefore cases with CNVs in feature-selected pathways, while “false positives” are controls with CNVs in feature-selected pathways. These definitions are used to compute the ROC curve below.

dat <- out$performance_denEnrichedNets
plot(0,0,type="n",xlim=c(0,100),ylim=c(0,100),
    las=1, xlab="False Positive Rate (%)", 
    ylab="True Positive Rate (%)",
    bty='n',cex.axis=1.5,cex.lab=1.3,
    main="ROC curve - Patients in label-enriched pathways")
points(dat$other_pct,dat$pred_pct,
      col="red",type="o",pch=16,cex=1.3,lwd=2)

plot of chunk unnamed-chunk-7

We can also compute the AUROC and AUPR from scratch.

tmp <- data.frame(  
    score=dat$score,
    tp=dat$pred_ol,fp=dat$other_ol,
    # tn: "-" that were correctly not called
    tn=dat$other_tot - dat$other_ol,
    # fn: "+" that were not called 
    fn=dat$pred_tot - dat$pred_ol) 

stats <- netDx::perfCalc(tmp)
tmp <- stats$stats
message(sprintf("PRAUC = %1.2f\n", stats$prauc))
## PRAUC = 0.67
message(sprintf("ROCAUC = %1.2f\n", stats$auc))
## ROCAUC = 0.69

Pathway scores are also added across the splits, for a total of 9 across the 3 splits (3 + 3 + 3).

# now get pathway score
print(head(out$cumulativeFeatScores))
##                                                                                                                             PATHWAY_NAME
## OVARIAN_TUMOR_DOMAIN_PROTEASES_cont.txt                                                                   OVARIAN_TUMOR_DOMAIN_PROTEASES
## NEUROTRANSMITTER_RECEPTORS_AND_POSTSYNAPTIC_SIGNAL_TRANSMISSION_cont.txt NEUROTRANSMITTER_RECEPTORS_AND_POSTSYNAPTIC_SIGNAL_TRANSMISSION
## HUNTINGTON_DISEASE_cont.txt                                                                                           HUNTINGTON_DISEASE
## L1CAM_INTERACTIONS_cont.txt                                                                                           L1CAM_INTERACTIONS
## NABA_ECM_AFFILIATED_cont.txt                                                                                         NABA_ECM_AFFILIATED
## TRANSCRIPTIONAL_REGULATION_OF_WHITE_ADIPOCYTE_DIFFERENTIATION_cont.txt     TRANSCRIPTIONAL_REGULATION_OF_WHITE_ADIPOCYTE_DIFFERENTIATION
##                                                                          SCORE
## OVARIAN_TUMOR_DOMAIN_PROTEASES_cont.txt                                      8
## NEUROTRANSMITTER_RECEPTORS_AND_POSTSYNAPTIC_SIGNAL_TRANSMISSION_cont.txt     7
## HUNTINGTON_DISEASE_cont.txt                                                  6
## L1CAM_INTERACTIONS_cont.txt                                                  6
## NABA_ECM_AFFILIATED_cont.txt                                                 6
## TRANSCRIPTIONAL_REGULATION_OF_WHITE_ADIPOCYTE_DIFFERENTIATION_cont.txt       6