R/PSEkNUCdi_RNA.R

Defines functions PSEkNUCdi_RNA

Documented in PSEkNUCdi_RNA

#' Pseudo k riboNucleotide Composition-Di(Parallel) (PSEkNUCdi_RNA)
#'
#' This function calculates the pseudo-k ribonucleotide composition(Di) (Parallel)
#' for each sequence.
#'
#' @details This function computes the pseudo ribonucleotide composition for each physicochemical property of di-ribonucleotides.
#' We have provided users with the ability to choose among the 22 properties in the di-ribonucleotide index database.
#'
#' @param seqs is a FASTA file containing ribonucleotide sequences. The sequences start
#' with '>'. Also, seqs could be a string vector. Each element of the vector is a ribonucleotide sequence.
#'
#'
#' @param selectedIdx is a vector of Ids or indices of the desired physicochemical properties of di-ribonucleotides.
#' Users can choose the desired indices by their ids or their names in the DI_RNA peoperties file.
#' The default value of this parameter is a vector with ("Rise (RNA)", "Roll (RNA)", "Shift (RNA)", "Slide (RNA)", "Tilt (RNA)","Twist (RNA)") ids.
#'
#' @param lambda is a tuning parameter. This integer value shows the maximum limit of spaces between di-ribonucleotide pairs. The Number of spaces
#' changes from 1 to lambda.
#'
#' @param w (weight) is a tuning parameter. It changes in the range of 0 to 1. The default value is 0.5.
#'
#' @param l This parameter keeps the value of l in lmer composition. The lmers form the first 4^l elements of the APkNCdi descriptor.
#'
#' @param ORF (Open Reading Frame) is a logical parameter. If it is set to true, ORF region of each sequence is considered instead of the original sequence (i.e., 3-frame).
#'
#' @param reverseORF is a logical parameter. It is enabled only if ORF is true.
#' If reverseORF is true, ORF region will be searched in the sequence and also in the reverse complement of the sequence (i.e., 6-frame).
#'
#' @param threshold is a number between (0 , 1]. In selectedIdx, indices with a correlation
#' higher than the threshold will be deleted. The default value is 1.
#'
#' @param label is an optional parameter. It is a vector whose length is equivalent to the number of sequences. It shows the class of
#' each entry (i.e., sequence).
#'
#'
#' @return a feature matrix such that the number of columns is 4^l+lambda and the number of rows is equal to the number of sequences.
#'
#'
#'
#' @export
#'
#' @examples
#'
#' fileLNC<-system.file("extdata/Carica_papaya101RNA.txt",package="ftrCOOL")
#' mat<-PSEkNUCdi_RNA(seqs=fileLNC,l=2,ORF=TRUE,threshold=0.8)


PSEkNUCdi_RNA<-function(seqs,selectedIdx=c("Rise (RNA)", "Roll (RNA)", "Shift (RNA)", "Slide (RNA)", "Tilt (RNA)","Twist (RNA)"),lambda=3,w = 0.05,l=2,ORF=FALSE,reverseORF=TRUE,threshold=1,label=c()){

  path.pack=system.file("extdata",package="ftrCOOL")

  if(length(seqs)==1&&file.exists(seqs)){
    seqs<-fa.read(seqs,alphabet="rna")
    seqs_Lab<-alphabetCheck(seqs,alphabet = "rna",label)

    seqs<-seqs_Lab[[1]]
    label<-seqs_Lab[[2]]
  }
  else if(is.vector(seqs)){
    seqs<-sapply(seqs,toupper)

    seqs_Lab<-alphabetCheck(seqs,alphabet = "rna",label)

    seqs<-seqs_Lab[[1]]
    label<-seqs_Lab[[2]]

  }
  else {
    stop("ERROR: Input sequence is not in the correct format. It should be a FASTA file or a string vector.")
  }
  flag=0
  if(ORF==TRUE){
    if(length(label)==length(seqs)){
      names(label)=names(seqs)
      flag=1
    }
    seqs=maxORF_RNA(seqs,reverse=reverseORF)
    if(flag==1)
      label=label[names(seqs)]
  }
  dict<-list("A"=1,"C"=2,"G"=3,"U"=4)
  numSeqs<-length(seqs)



  #aaIdxAD<-paste0(path.pack,"/standardizeNUCidx2.csv")
  aaIdxAD<-paste0(path.pack,"/DI_RNA.csv")
  aaIdx<-read.csv(aaIdxAD)

  row.names(aaIdx)<-aaIdx[,1]
  aaIdx<-aaIdx[selectedIdx,-1]
  aaIdx<-as.matrix(aaIdx)
  aaIdx<-type.convert(aaIdx)


  ###deleteing high correlate features
  if(threshold!=1){

    aaIdx<-t(aaIdx)
    corr<-cor(aaIdx)
    corr2<-corr^2
    tmp<-corr2
    tmp[upper.tri(tmp)]<-0
    for(i in 1:ncol(tmp)){
      tmp[i,i]=0
    }
    aaIdx<- aaIdx[,!apply(tmp,2,function(x) any(x > threshold))]
    aaIdx<- t(aaIdx)
  }##normalize beween [0,1]
  minFea<-apply(aaIdx, 1, min)
  maxFea<-apply(aaIdx, 1, max)
  aaIdx<-(aaIdx-minFea)/(maxFea-minFea)

  numFea<-nrow(aaIdx)

  start<-4^l+1
  end<-4^l+lambda

  featureMatrix<-matrix(0,nrow = numSeqs,ncol = ((4^l)+lambda))
  for(n in 1:numSeqs){
    seq<-seqs[n]
    N<-nchar(seq)
    if (lambda>N || lambda<=0){
      stop("Error: lambda should be between [1,N]. N is the minimum of sequence lengths")
    }


    small_theta<-vector(mode = "numeric", length = lambda)
    chars<-unlist(strsplit(seq,NULL))
    temp1<-chars[1:(N-1)]
    temp2<-chars[2:N]
    dimers<-paste0(temp1,temp2)
    lenDimer=N-1

    for(k in 1:lambda){
      vecti=1:(lenDimer-k)
      vectj=vecti+k
      tempMat<-matrix(0,nrow = numFea,ncol = (lenDimer-k))




      for(m in 1:numFea){

        tempMat[m,]=(aaIdx[m,dimers[vecti]]-aaIdx[m,dimers[vectj]])^2
      }

      bigThetaVect<-apply(tempMat, 2, sum)
      sumbigThetas<-sum(bigThetaVect)
      small_theta[k]<-(1/(N-k))*sumbigThetas
    }

    sum_small_th<-sum(small_theta)
    NUC<-kNUComposition_RNA(seq,rng=l,normalized = FALSE,upto = FALSE)
    featureVector <- vector(mode = "numeric", length = (length(NUC)+lambda))
    featureVector[1:length(NUC)]<-(NUC)/(N+w*(sum_small_th))

    featureVector[start:end]<- (w*small_theta)/(N+w*(sum_small_th))


    featureMatrix[n,]=featureVector


  }
  namNUC<-nameKmer(k=l,type = "rna")
  temp=1:lambda
  colnam=c(namNUC,paste("lambda",temp,sep=""))
  colnames(featureMatrix)=colnam
  if(length(label)==numSeqs){
    featureMatrix<-as.data.frame(featureMatrix)
    featureMatrix<-cbind(featureMatrix,label)
  }
  row.names(featureMatrix)<-names(seqs)
  return(featureMatrix)

}

Try the ftrCOOL package in your browser

Any scripts or data that you put into this service are public.

ftrCOOL documentation built on Nov. 30, 2021, 1:07 a.m.