.packageName <- "affypdnn"
expressopdnn <-  function(afbatch,
                          ## --
                          bg.correct=FALSE,
                          bgcorrect.method = NULL,
                          bgcorrect.param = list(),
                          ## --
                          normalize = FALSE,
                          normalize.method = NULL,
                          normalize.param=list(),
                          ## --
                          pmcorrect.method = c("pdnn", "pdnnpredict"),
                          ##pmcorrect.param = list(),
                          ## --
                          params.chiptype = NULL,
                          summary.subset = NULL,
                          ## ---
                          eset.normalize = TRUE,
                          scale.to = 500,
                          verbose = TRUE,
                     widget = FALSE
                     ) {

  if (widget) {
    require(tkWidgets) || stop("library tkWidgets could not be found !")
    warning("Option 'widget' not yet implemented ! (reverting the option to 'FALSE')")
    widget <- FALSE
  }

  nchips <- length(afbatch)
  
###background stuff must be added before normalization!
  
  ## -- background correction method
  if (bg.correct && is.null(bgcorrect.method)) {
    if (widget) {
      bgcorrect.method <- pwidget.selector(bgcorrect.methods,
                                    title = "Method for background correction")
                                    ##,choices.help = NULL)
      #bg.method <- paste("bg.correct", bg.method, sep=".")
    } else {
      stop("bg.method missing")
    }
  }

  ## -- normalize.method
  if ((normalize) & (is.null(normalize.method))) {
    if (widget) {
      ##DEBUG: move what's below to 'pwidget.selector'
      n.methods <- normalize.methods(afbatch)
      normalize.method <- pwidget.selector(n.methods,
                                           title = "Method for normalization")
      ##choices.help = paste("normalize", n.methods,  sep="."))##took this out cause netscape popping up was annoying me
                                                                                
      rm(n.methods)
    } else {
      stop(paste("normalization method missing. Try one of:",
                 normalize.methods(afbatch), sep="\n"))
    }
  }

    ## -- pm correction method
  if (is.null(pmcorrect.method)) {
    if (widget) {
      pmcorrect.method <- pwidget.selector(c("pdnn", "pdnnpredict"),
                                    title = "Method for PM correction")
      ##,choices.help = NULL)
      #bg.method <- paste("bg.correct", bg.method, sep=".")
    } else {
      stop("pmcorrect.method missing")
    }
  } else {
    pmcorrect.method <- match.arg(pmcorrect.method)
  }
  
  ## -- summary of what will be done
  if (verbose) {
    if (bg.correct){
      cat("background correction:", bgcorrect.method, "\n")
    }
    if (normalize) {
      cat("normalization:", normalize.method, "\n")
    }
    cat("PM/MM correction :", "pdnn", "\n")
    cat("expression values:", "pdnn", "\n")
  }
  
    ## -- background correct (if needed)
  if (bg.correct) {
    
    if (verbose)
      cat("background correcting...")
    
    afbatch <- do.call("bg.correct", c(alist(afbatch, method=bgcorrect.method),
                                       bgcorrect.param))
    
    if (verbose)
      cat("done.\n")
  }


  ## -- normalize (if wished)
  if (normalize) {
                                                                                
    if (verbose)
      cat("normalizing...")
    
    afbatch <- do.call("normalize",
                       c(alist(afbatch, normalize.method), normalize.param))
    
    if (verbose)
      cat("done.\n")
  }

  ## chip-type specific parameters
  if (is.null(params.chiptype)) {
    ## try to get it from the pack
    namebase <- cleancdfname(afbatch@cdfName) 
    dataname <- paste(substr(namebase, 1,  nchar(namebase) - 3), ".pdnn.params", sep="")
    if(! dataname %in% do.call("data", list(package="affypdnn"))$results[, 3])
      stop("params.chiptype missing !")
    do.call("data", list(dataname, package="affypdnn"))
    assign("params.chiptype", get(dataname))
  }
  
  params <- find.params.pdnn(afbatch, params.chiptype)
  
  eset <- computeExprSet(afbatch, pmcorrect.method=pmcorrect.method,
                         summary.method="pdnn",
                         ids=summary.subset,
                         pmcorrect.param = list(params, params.chiptype=params.chiptype, callingFromExpresso=TRUE),
                         summary.param = list(params))

  if (eset.normalize)
    eset <- pdnn.scalevalue.exprSet(eset, scale.to=scale.to)
  
  return(eset)
  
}
find.params.pdnn <- function(abatch, params.chiptype, optim.method="BFGS", verbose=TRUE, give.warnings=TRUE) {

  if (verbose)
    cat("initializing data structure...")

  names.abatch <- sort(geneNames(abatch))
  n.probes <- rep(0, length(names.abatch))
  n.cel <- length(sampleNames(abatch))
  ##names.pdnnparams <- sort(ls(params.chiptype$params.gene))
  names.i <- rep(as.integer(NA), length(names.abatch))
  lambda <- vector("list", length=length(names.abatch))
  Fs <- rep(as.numeric(NA), length=n.cel)
  S.lambda.E <- vector("numeric", length=length(names.abatch)) ## \sum \frac{\lambda_{ij}} {1+exp(E_{ij})}
  
  ##K1top <- vector("numeric", length=length(names.abatch))
  ##K2top <- vector("numeric", length=length(names.abatch))
  ##K3top <- vector("numeric", length=length(names.abatch))
  ##prealpha <- alpha <- beta <- rep(as.numeric(NA), length(names.abatch))
  Bs <- vector("numeric", length=n.cel)
  Ns <- vector("numeric", length=n.cel)

  ## init data structures
  for (gene.i in seq(along=names.abatch)) {
    gene <- names.abatch[gene.i]
    if(! exists(gene, envir=params.chiptype$params.gene)) {
      if (give.warnings) {
        warning(paste("The gene", gene, "could not be found in the parameter files (possible entanglement with cdfenvs) !"))
      }
      n.probes[gene.i] <- 0
      names.abatch[gene.i] <- NA
    } else {
      names.i[gene.i] <- get(gene, envir=params.chiptype$params.gene)
      n.probes[gene.i] <- length(params.chiptype$gene.Sn[[names.i[gene.i]]])
      lambda[[gene.i]] <- matrix(as.numeric(NA), n.probes[gene.i], n.cel)
    }
    S.lambda.E[gene.i] <- NA
  }

  names.abatch.nonNA <- which(! is.na(names.abatch))
  i.pm <- indexProbes(abatch, which="pm", genenames=names.abatch)
  all.i.pm <- unlist(i.pm[names.abatch.nonNA])

  if (verbose)
    cat("done.\n")
  
  ## loop across chips
  for (cel.i in seq(along=sampleNames(abatch))) {
    
    if (verbose)
      cat("dealing with CEL", cel.i, ":\n")

    if (verbose)
      cat("  step 1...")
    
    for (gene.i in names.abatch.nonNA) {

      gene <- names.abatch[gene.i]
      
      gene.params.i <- names.i[gene.i]
      gene.Sg <- params.chiptype$gene.Sg[[gene.params.i]]
      gene.Sn <- params.chiptype$gene.Sn[[gene.params.i]]
      gene.xy <- params.chiptype$gene.xy[[gene.params.i]]
      gene.o <- match(i.pm[[gene.i]], xy2indices(gene.xy[,1], gene.xy[,2], abatch=abatch))
      probe.intensities <- intensity(abatch)[i.pm[[gene.i]], cel.i, drop=FALSE]
      
      lambda[[gene.i]][, cel.i] <- sqrt(probe.intensities * gene.Sg[gene.o])
      S.lambda.E[gene.i] <- sum( sqrt(probe.intensities / gene.Sg[gene.o]) )
      ##this.lambda <- lambda[[gene.i]][, cel.i]
      ##lambda[[gene]] <- sqrt(intensity(abatch)[i.pm, ] * get(gene, envir=params.chiptype$Sg))
       ## prealpha[gene.i] <- sum(1 / (gene.Sg[gene.o] * this.lambda))
      ##K1top[gene.i] <- sum(probe.intensities / this.lambda)
      ##K2top[gene.i] <- sum(1 / this.lambda)
      ##K3top[gene.i] <- sum(1/(gene.Sn[gene.o] * this.lambda))
    }
    
    if (verbose)
      cat("done.\n")

    ##K1top.rep <- rep(K1top, n.probes)
    ##K2top.rep <- rep(K2top, n.probes)
    ##K3top.rep <- rep(K3top, n.probes)
    
    ##all.gene.params <- multiget(names.abatch, envir=params.chiptype$params.gene)
    ##alpha <- rep(prealpha, n.probes) * unlist(lapply(all.gene.params, function(x) x$Sg))
    ##KI1 <- (K1top.rep / alpha) / intensity(abatch)[all.i.pm, cel.i, drop=FALSE]
    ##KI2 <- (K2top.rep / alpha - 1) / intensity(abatch)[all.i.pm, cel.i, drop=FALSE]
    ##KI3 <- (K3top.rep / alpha - (1 / unlist(lapply(all.gene.params, function(x) x$Sn)))) /
    ##  intensity(abatch)[all.i.pm, cel.i, drop=FALSE]
    ##rm(all.gene.params, K1top.rep, K2top.rep, K3top.rep)
    
    if (verbose)
      cat("  step 2...")
    ##FIXME: starting values always the same ?
    B<-50
    N<-3000

    lambda.cel <- unlist(lapply(names.abatch.nonNA, function(x) lambda[[x]][,cel.i]))
    Betaij.left <- sum(1 / unlist(params.chiptype$gene.Sg[seq(along=names.abatch.nonNA)]) / lambda.cel)
    Betaij <- lapply(seq(along=names.abatch.nonNA), function(i, x, y) y[[i]] * x, Betaij.left, params.chiptype$gene.Sg)
    Gij.top <- sum(intensity(abatch)[all.i.pm, cel.i] / lambda.cel)
    Gij <- lapply(seq(along=names.abatch.nonNA), function(i, x, y) x / y[[i]], Gij.top, Betaij)
    Hij.top <- sum(1 / lambda.cel)
    Hij <- lapply(seq(along=names.abatch.nonNA), function(i, x, y) x / y[[i]] - 1, Hij.top, Betaij)
    Kij.lefttop <- sum(1 / lambda.cel / unlist(params.chiptype$gene.Sn[seq(along=names.abatch.nonNA)]))
    Kij <- lapply(seq(along=names.abatch.nonNA), function(i, x, y, z) x / y[[i]] + 1 / z[[i]], Kij.lefttop, Betaij, params.chiptype$gene.Sn)
    intensity.forfit <- intensity(abatch)[unlist(i.pm[names.abatch.nonNA]), cel.i]
    intensity.forfit.positive <- intensity.forfit > 0
    Gij.forfit <- unlist(Gij[names.i[names.abatch.nonNA]])
    Hij.forfit <- unlist(Hij[names.i[names.abatch.nonNA]])
    Kij.forfit <- unlist(Kij[names.i[names.abatch.nonNA]])
    
    ##fit.f<-function(par, Gij, Hij, Kij, intensities.matrix, i.pm){
    ##fit.f<-function(par){
    fit.f<-function(par, Gij.forfit, Hij.forfit, Kij.forfit) {
      B <- par[1]
      N <- par[2]
      sum.sqrt <- 0
      n.probes <- 0
      Ihat <-  Gij.forfit - B * Hij.forfit - N * Kij.forfit
      good <- Ihat > 0 & intensity.forfit.positive
      sum.sqr <- sum((log(Ihat[good] / intensity.forfit[good]))^2)
#       for (gene.i in names.abatch.nonNA) {
#         ##cat(gene.i, "\n")
#         i <- names.i[gene.i]
#         probe.intensities <- intensity.matrix[i.pm[[gene.i]], cel.i]
#         I.hat <- Gij[[i]] - B * Hij[[i]] - N * Kij[[i]]
#         good <- I.hat > 0 & probe.intensities > 0
#         n.probes <- n.probes + sum(good)
#         sum.sqr <- sum.sqrt + sum((log(I.hat[good]) - log(probe.intensities[good]))^2)
#       }
      ##cat(sum.sqr / sum(good), "\n")
      return(sum.sqr / sum(good))
      ##return( sum(log(Y[Y >= 0])**2) / length(Y[Y >= 0]))
    }

    if (optim.method == FALSE) {
      ## use steepest descent
      i<-succes<-failed<-0

      DeltaN<-500;DeltaB<-10
      Fit <- +Inf
      Bdirection<-1
      Ndirection<-1
      improvement<-2
      
      while (improvement>1.00000001 || DeltaN>0.5 || DeltaB>0.5){
        DB<-DeltaB*runif(1,0.85,1)
        Bnew<-B+Bdirection*DB
        DN<-DeltaN*runif(1,0.85,1)
        Nnew<-N+Ndirection*DN
        Fnew <- fit.f(c(Bnew,Nnew), Gij.forfit=Gij.forfit, Hij.forfit=Hij.forfit, Kij.forfit=Kij.forfit)
        i<-i+1
        if(Fnew<Fit){# Succes
                                        #cat("PDNN: Iteration:",i, " N':", N, " B:", B, " F:", Fnew, "\n")#(Ndirection*DeltaN)/(Bdirection*DeltaB)
          improvement <- Fit/Fnew
          succes<-succes+1
          failed<-0
          DeltaB<-DB; DeltaN<-DN
          Fit<-Fnew; B<-Bnew; N<-Nnew
          if(succes>1){
            increase<-runif(1,1,1.5)
            DeltaB<-DeltaB*increase
            DeltaN<-DeltaN*increase
          }
          
        }else{ # bad guess
          failed<-failed+1
          succes<-0
          if(failed>1){
            DeltaB<-DeltaB*runif(1,0.85**failed,1)
            DeltaN<-DeltaN*runif(1,0.9**failed,1)
            if(runif(1,0,1)>0.4999){Bdirection<-1}else{Bdirection<--1}
            if(runif(1,0,1)>0.4999){Ndirection<-1}else{Ndirection<--1}
          }
        }
      }
      Bs[cel.i] <- B
      Ns[cel.i] <- N
    } else {
      ## use R's 'optim'
      
      ##par.f <- optim(c(B, N), fit.f, method=optim.method, lower=-Inf, upper=-Inf, control=list(), hessian=FALSE, Betaij, Gij, Hij, Kij)
      par.f <- optim(c(B, N), fit.f, method=optim.method, Gij=Gij.forfit, Hij=Hij.forfit, Kij=Kij.forfit)
      
      Bs[cel.i] <- par.f$par[1]
      Ns[cel.i] <- par.f$par[2]
    }
    if (verbose)
      cat("done.\n")

    Fs[cel.i] <- fit.f(c(Bs[cel.i], Ns[cel.i]), Gij.forfit, Hij.forfit, Kij.forfit)
    
  }


  
  return(list(lambda=lambda, Bs=Bs, Ns=Ns, Fs=Fs, names.abatch=names.abatch, names.i=names.i))
}
matplotProbesPDNN <- function(x, type="l", ...) {
  matplot(x, type=type, ...)  
  ok <- attr(x, "ok")
  for (i in seq(1, ncol(x), length=ncol(x))) {
    points(x[, i], pch=c(1,3)[as.integer(ok[, i])+1])
  }
}

pmcorrect.pdnnpredict <- function(object, params, gene=NULL, gene.i=NULL, params.chiptype=NULL, outlierlim=3, callingFromExpresso=FALSE) {  

  new.int <- pmcorrect.pdnn(object, params, gene=gene, gene.i=gene.i, params.chiptype=params.chiptype, outlierlim=outlierlim, callingFromExpresso=callingFromExpresso)
  ## do the prediction bit
  
  Bs <- params$Bs
  Ns <- params$Ns
  Sn.gene <- attr(new.int, "Sn.gene")

 
  new.int <- sweep(new.int + 1/outer(Sn.gene, Ns, "/"), 2, Bs, "+")
  
  return(new.int)
}

pmcorrect.pdnn <- function(object, params, gene=NULL, gene.i=NULL, params.chiptype=NULL, outlierlim=3, callingFromExpresso=FALSE) {  

  probes <- object@pm
  if (is.null(params.chiptype)) {
    stop("params.chiptype must be specified.")
  }

  lambda <- params$lambda
  Bs <- params$Bs
  Ns <- params$Ns
  Fs <- params$Fs

  if (is.null(gene.i)) {
    if (callingFromExpresso) {
      gene <- get("id", envir=parent.frame(3)) ## dynamic lexical scoping... (not static)
    } else {
      if (is.null(gene))
        ##stop("gene.i must be specified.")
        gene <- object@id
    }
  }

  gene.i <- get(gene, envir=params.chiptype$params.gene)
  Sg.gene <- params.chiptype$gene.Sg[[gene.i]]
  Sn.gene <- params.chiptype$gene.Sn[[gene.i]]
  
  Ntop <- sweep(probes, 2, Bs , "-") - sapply(Ns, "/", Sn.gene)
  
  new.int <- probes
  ok <- matrix(as.logical(NA), nr=nrow(probes), nc=ncol(probes))
  
  for (cel.i in seq(1, ncol(probes), length=ncol(probes))) {
    new.int[, cel.i] <- sum((Ntop / lambda[[gene.i]][, cel.i]) / sum(1 / Sg.gene / lambda[[gene.i]][, cel.i])) /
      Sg.gene #+ Ns[cel.i] / Sn.gene + Bs[cel.i]
    
    ok[, cel.i] <- Ntop[, cel.i] >= 0 & probes[, cel.i] / new.int[, cel.i] > 0 & log(probes[, cel.i] / new.int[, cel.i]) < outlierlim * sqrt(Fs[cel.i]) & !is.na(log(probes[, cel.i] / new.int[, cel.i]))
  }
  attr(new.int, "ok") <- ok
  attr(new.int, "Ntop") <- Ntop
  attr(new.int, "gene.i") <- gene.i
  attr(new.int, "Sg.gene") <- Sg.gene
  attr(new.int, "Sn.gene") <- Sn.gene

  return(new.int)
}


generateExprVal.method.pdnn <- function(probes, params) {  
  
#   if (is.null(params.chiptype)) {
#     stop("params.chiptype must be specified.")
#   }

  if (is.null(attr(probes, "ok")) || is.null(attr(probes, "Ntop")) || is.null(attr(probes, "gene.i")))
    stop("The summary method 'pdnn' can only be used with the pmcorrect method 'pdnn' !")

  lambda <- params$lambda
#   Bs <- params$Bs
#   Ns <- params$Ns
#   Fs <- params$Fs

  
#   if (is.null(gene.i)) {
#     if (callingFromExpresso)
#       gene <- get("id", envir=parent.frame(3)) ## dynamic lexical scoping... (not static)
#     else
#       stop("gene.i must be specified.")
#   }    
#   params.gene <- get(gene, envir=params.chiptype$params.gene)
#   gene.i <- params.gene$gene.i
#   Sg.gene <- params.gene$Sg
#   Sn.gene <- params.gene$Sn

  ##gene <- names.abatch[gene.i]
  ##names(Nj) <- names(Iobs)
  ##options(warn=-1)

  ##FIXME: useless ?
  ##i.pm <- indexProbes(abatch, gene, which="pm")
  ##i.mm <- indexProbes(abatch, gene, which="mm")

  ##cat(str(probes))
  ##cat(str(Bs+Ns))
  ##cat(str(get(gene, envir=params.chiptype$Sn)))
  
#   Ntop <- sweep(probes, 2, Bs , "-") - sapply(Ns, "/", params.gene$Sn)
  
#   new.int <- probes
#   ok <- matrix(as.logical(NA), nr=nrow(probes), nc=ncol(probes))
  
#   for (cel.i in seq(1, ncol(probes), length=ncol(probes))) {
#     new.int[, cel.i] <- sum((Ntop / lambda[[gene.i]][, cel.i]) / sum(1 / Sg.gene / lambda[[gene.i]][, cel.i])) /
#       Sg.gene + Ns[cel.i] / Sn.gene + Bs[cel.i]
#     ok[, cel.i] <- Ntop[, cel.i] >= 0 & probes[, cel.i] / new.int[, cel.i] > 0 & log(probes[, cel.i] / new.int[, cel.i]) < outlierlim * sqrt(Fs[cel.i]) & !is.na(log(probes[, cel.i] / new.int[, cel.i]))
#   }

  Ntop <- attr(probes, "Ntop")
  Sg.gene <- attr(probes, "Sg.gene")
  gene.i <- attr(probes, "gene.i")
  ok <- attr(probes, "ok")
  
  expr.val <- rep(as.numeric(NA), ncol(probes))
  expr.se <- rep(as.numeric(NA), ncol(probes))
  
  for (cel.i in seq(1, ncol(probes), length=ncol(probes)))
    expr.val[cel.i] <- sum((Ntop / lambda[[gene.i]][, cel.i])[ok[, cel.i]]) / sum((1/Sg.gene / lambda[[gene.i]][, cel.i])[ok[, cel.i]])
       
  return(list(exprs=expr.val, se.exprs=expr.se))
}

pdnn.scalevalue.exprSet <- function(eset, scale.to=500) {
  m <- exprs(eset)
  m.mean <- apply(m, 1, mean, na.rm=TRUE)
  mm <- sweep(m, 2, scale.to/m.mean, "*")
  exprs(eset) <- mm
  return(eset)
}
pdnn.params.chiptype <- function(energy.param.file, probes.file = NULL, probes.pack = NULL, probes.data.frame = NULL,
                                 seq.name="sequence", x.name="x", y.name="y", affyid.name="Probe.Set.Name",
                                 verbose=TRUE) {


  if (! file.exists(energy.param.file))
    stop(paste("Cannot find the energy.param.file '", energy.param.file, sep=""))
  
  if (sum(c(is.null(probes.file), is.null(probes.pack), is.null(probes.data.frame))) != 2) {
    stop("Specify one 'probe.file' _or_ one 'probe.pack' _or_ 'probe.data.frame'")
  }

#   if (is.null(cdfName)) {
#     stop("'cdfName', a name to find the corresponding cdfenv is missing !")
#   }

#   ## Hack for version 1.2.x of 'affy'
#   a <- new("AffyBatch", cdfName=cdfName)
#   cdfenv <- getCdfInfo(a)
#   rm(a)
  
  ## FIXME
  if (!is.null(probes.file)) {
    probe.tab <- read.table(probes.file, sep="\t",header=TRUE, nrows=2)##[c(1,2,3,5)]
  }
  
  if (!is.null(probes.pack)) {
    do.call("library", list(probes.pack))
    probe.tab <- get(probes.pack, envir=as.environment(paste("package:", probes.pack, sep="")))
  }

  if (!is.null(probes.data.frame))
    probe.tab <- probes.data.frame
  ##
  
  i.seq <- match(seq.name, names(probe.tab))
  i.x <- match(x.name, names(probe.tab))
  i.y <- match(y.name, names(probe.tab))
  i.affyid <- match(affyid.name, names(probe.tab))

  if (any(is.na(c(i.seq, i.x, i.y, i.affyid)))) {
    m <- paste("You asked :", paste(" affyid.name:", affyid.name), paste(" x.name  :", x.name),
               paste(" y.name  :", y.name), paste(" seq.name :", seq.name),
               paste("While the names in the data.frame are:"),
               paste(paste("",names(probe.tab)), collapse="\n"), sep="\n")
    stop(m)
  }

  if (! is.null(probes.file)) {
    if (verbose)
      cat("Reading the probes file...")
    probe.tab <- read.table(probes.file, sep="\t",header=TRUE)
    if (verbose)
      cat("done.\n")
  }
  
  probe.x <- probe.tab[[i.x]] + 1
  probe.y <- probe.tab[[i.y]] + 1
  affy.id <- probe.tab[[i.affyid]]
  probe.seq <- tolower(as.character(probe.tab[[i.seq]]))
  
  if(verbose)
    cat("Calculating chip type specific parameters, (may take some time)...\n")

  ## Reading probe sequence, energy and position weight files
  
  ## FIXME (automagic to do) ?
  
  ep <- read.table(energy.param.file, nrows=80, header=TRUE)

  Wg <- as.vector(ep[33:56, 2]) ## weights (specific)
  
  Wn <- as.vector(ep[57:80, 2]) ## weights (non-specific)

#    params.chiptype <- list(Eg = new.env(hash=TRUE),
#                            Wg = Wg,
#                            En= new.env(hash=TRUE),
#                            Wn = Wn,
#                            Sg = new.env(hash=TRUE),
#                            Sn = new.env(hash=TRUE),
#                            gene2i = new.env(hash=TRUE))
  
  params.chiptype <- list(Eg = new.env(hash=TRUE),
                          Wg = Wg,
                          En= new.env(hash=TRUE),
                          Wn = Wn,
                          gene.Sn = NA,
                          gene.Sg = NA,
                          gene.xy = NA,
                          gene.name = NA, ## only here for cross-checks
                          params.gene = new.env(hash=TRUE))
  
  
  Eg <- multiassign(as.list(ep[[1]][1:16]), as.list(as.vector(ep[1:16, 2])), envir=params.chiptype$Eg) ## energy (specific)
  
  ##Wg <- multiassign(rownames(ep), as.vector(ep[33:56, 2]), envir=env.list$Wg) ## weights (specific)
  
  En <- multiassign(as.list(ep[[1]][17:32]), as.list(as.vector(ep[17:32, 2])), envir=params.chiptype$En) ## energy (non-specific)
  
  
  ##Wg <- params.chiptype$Wg
  ##Wn <- params.chiptype$Wn
  En <- params.chiptype$En
  Eg <- params.chiptype$Eg
  
  Sf <- function(oligo){
    ## oligo: vector of probe sequences
    ## calculating (1/exp(E)) and (1/exp(E*)

    ENv <-EGv <-vector(length = length(oligo))

    increment <- as.integer(1)
    ## loop across the oligos
    for (g in 1:length(oligo)){
      
      EG <- EN <- 0

      ## walk along the sequence
      for (k in seq(1, nchar(oligo[g])-1)) {

        ## FIXME: build hashtables for Eg and En
        ##di.nucl <- substr(oligo[g], k, k+1)
        di.nucl <- .Internal(substr(oligo[g], k, k+increment))
        ## FIXME cast necessary ?
        ##EG <- as.numeric(EG + Wg[k] * get(di.nucl, envir = Eg))
        ##EN <- as.numeric(EN + Wn[k] * get(di.nucl, envir = En))
        EG <- EG + Wg[k] * get(di.nucl, envir = Eg)
        EN <- EN + Wn[k] * get(di.nucl, envir = En)
        
      }
      
      EGv[g] <- 1 + exp(EG)
      ENv[g] <- 1 + exp(EN)
    }
    
    return(cbind(EGv, ENv))
  }

  
  ## FIXME (speedup w/ affy cdfenvs ?)
  ##S <- tapply(as.vector(probe.tab[[i.seq]]), affy.id, Sf)
  S.index <- tapply(seq(1, nrow(probe.tab), length=nrow(probe.tab)), affy.id, I, simplify=FALSE)
  
  params.chiptype$gene.Sn <- vector("list", length=length(S.index))
  params.chiptype$gene.Sg <- vector("list", length=length(S.index))
  params.chiptype$gene.xy <- vector("list", length=length(S.index))
  params.chiptype$gene.name <- vector("character", length=length(S.index))
  
  params.gene <- params.chiptype$params.gene

  if (verbose) {
    pbt <- new("ProgressBarText", length(S.index), barsteps = as.integer(20))
    open(pbt)
  }
  
  ## FIXME:
  for (gene.i in seq(along=S.index)) {

    if (verbose)
      update(pbt)
    
    gene <- names(S.index)[gene.i]
    gene.S.index <- S.index[[gene.i]]
    sequences.in.ppset <- probe.seq[gene.S.index]
    S <- Sf(sequences.in.ppset)
    ## FIXME store the Xs and Ys for now (better scheme when cdfenvs in classes)
    
    ##assign(gene, S[[gene.i]][, 1], envir=Sg) #Sg is now (1+exp(E))
    ##assign(gene, S[[gene.i]][, 2], envir=Sn) #Sn is now (1+exp(E*))
    ##assign(gene, gene.i, envir=params.chiptype$gene2i)
    assign(gene, gene.i, envir=params.gene)
    params.chiptype$gene.Sg[[gene.i]] <- S[, 1]
    params.chiptype$gene.Sn[[gene.i]] <- S[, 2]
    params.chiptype$gene.xy[[gene.i]] <- cbind(probe.x[gene.S.index], probe.y[gene.S.index])
    params.chiptype$gene.name[[gene.i]] <- gene
  }
  
  if (verbose)
    close(pbt)
  
  return(params.chiptype)
  
}

transform.ProbeSet <- function(x, fun=I, ...) {
  x@pm <- fun(x@pm, ...)
  x@mm <- fun(x@mm, ...)
  return(x)
}
.First.lib <- function(libname, pkgname, where) {
  
  require(affy)

  ## register the new summary method
  where <- match(paste("package:", pkgname, sep=""), search())
  cat("registering new summary method 'pdnn'.\n")
  assign("express.summary.stat.methods", c(express.summary.stat.methods, "pdnn"), envir=as.environment("package:affy"))
  cat("registering new pmcorrect method 'pdnn' and 'pdnnpredict'.\n")
  assign("pmcorrect.methods", c(pmcorrect.methods, "pdnn", "pdnnpredict"), envir=as.environment("package:affy"))
  ##assign("express.summary.stat.methods", c(express.summary.stat.methods, "pdnn"), envir=where)
  
}
