.packageName <- "Bhat"
"btrf" <-
function(xt,xl,xu) {

####  back transformation
####  this assumes logit transformations of the parameters
####  bounded from below by xl and from above by xu

  rho <- exp(xt)
  return((rho * xu + xl)/(1.+rho))
}
"dfp" <-
function (x, f, tol=1e-5, nfcn = 0) 
{
	    #     Function Minimization for R. 
	    #     This function is part of the Bhat exploration tool and is
	    #     based on Fletcher's "Switching Method" (Comp.J. 13,317 (1970)).

	    #     This program is free software; you can redistribute it and/or modify
	    #     it under the terms of the GNU General Public License as published by
	    #     the Free Software Foundation; either version 2 of the License, or
	    #     (at your option) any later version.
	    #     This program is distributed in the hope that it will be useful,
	    #     but WITHOUT ANY WARRANTY; without even the implied warranty of
	    #     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
	    #     GNU General Public License for more details.

            #     E Georg Luebeck (gluebeck@fhcrc.org)

            xt.inf <- 16
	    slamin <- 0.2
	    slamax <- 3
	    tlamin <- 0.05
	    tlamax <- 6
	    iter <- 0
            status <- 1
	    if (!is.list(x)) {
	    	    cat("x is not a list! see help file", "\n")
	    	    return()
	    }
	    names(x)[1] <- "label"
	    names(x)[2] <- "est"
	    names(x)[3] <- "low"
	    names(x)[4] <- "upp"
	    npar <- length(x$est)
	    ####  objects:
	    if (npar <= 0) {
	    	    warning("no. of parameters < 1")
	    	    exit
	    }
	    g <- numeric(npar)
	    gs <- numeric(npar)
	    g2 <- numeric(npar)
	    xt <- numeric(npar)
	    xxs <- numeric(npar)
	    dirin <- numeric(npar)
	    delgam <- numeric(npar)
	    vg <- numeric(npar)
	    ####  first call 
	    cat(date(), "\n","\n")
	    xt <- ftrf(x$est, x$low, x$upp)
	    fmin <- f(x$est)
	    nfcn <- nfcn + 1
            cat('starting at','\n')
	    cat(format(nfcn), "  fmin: ", fmin, "   ", format(x$est), "\n")
	    isw2 <- 0
	    nretry <- 0
	    nfcnmx <- 10000
	    vtest <- 0.001
	    apsi <- 0.05
	    ####  or .01 for better convergence?
	    up <- 1
	    ####  memorize current no. of function calls
	    npfn <- nfcn
	    rho2 <- 10 * apsi
	    rostop <- tol * apsi
	    trace <- 1
	    fs <- fmin
	    cat("\n")
	    # cat("rostop: ", rostop, "\n", "apsi:   ", apsi, "\n", "vtest: ", vtest, "\n")
	    dfp.loop <- 0
	    # while()
	    while (dfp.loop < 10) {
	    	    #### 1, 2 
	    	    #### ?? cat("COVARIANCE MATRIX NOT POSITIVE-DEFINITE","\n")
	    	    #### define step sizes dirin
                    d <- dqstep(list(label=x$label,est=btrf(xt, x$low, x$upp),low=x$low,upp=x$upp),f,sens=.01)
	    	    if (isw2 >= 1)  d <- 0.02 * sqrt(abs(diag(v)) * up)
	    	    dirin <- d

	    	    ##### obtain gradients, second derivatives, search for pos. curvature #########
	    	    ntry <- 0
	    	    negg2 <- 1
	    	    # loop 10

	    	    while (negg2 >= 1 & ntry <= 6) {
	    	    	    negg2 <- 0
	    	    	    #for(id in 1:npar)  loop 10
	    	    	    for (i in 1:npar) {
	    	    	    	    d <- dirin[i]
	    	    	    	    xtf <- xt[i]
	    	    	    	    xt[i] <- xtf + d
	    	    	    	    fs1 <- f(btrf(xt, x$low, x$upp))
	    	    	    	    nfcn <- nfcn + 1
	    	    	    	    xt[i] <- xtf - d
	    	    	    	    fs2 <- f(btrf(xt, x$low, x$upp))
	    	    	    	    nfcn <- nfcn + 1
	    	    	    	    xt[i] <- xtf
	    	    	    	    gs[i] <- (fs1 - fs2)/(2 * d)
	    	    	    	    g2[i] <- (fs1 + fs2 - 2 * fmin)/d^2
	    	    	    	    #if (g2[i]  <=  0.) 
	    	    	    	    if (g2[i] <= 0) {
	    	    	    	     ####  search if g2  <=  0. . .
	    	    	    	     cat("covariance matrix is not positive-definite or nearly singular", 
	    	    	    	      "\n")
	    	    	    	     negg2 <- negg2 + 1
	    	    	    	     ntry <- ntry + 1
	    	    	    	     d <- 50 * abs(dirin[i])
	    	    	    	     xbeg <- xtf
	    	    	    	     if (gs[i] < 0) {
	    	    	    	      dirin[i] <- -dirin[i]
	    	    	    	     }
	    	    	    	     kg <- 0
	    	    	    	     nf <- 0
	    	    	    	     ns <- 0
	    	    	    	     # while(ns < 10 ...)
	    	    	    	     while (ns < 10 & nf < 10) {
	    	    	    	      xt[i] <- xtf + d
	    	    	    	      f0 <- f(btrf(xt, x$low, x$upp))
	    	    	    	      nfcn <- nfcn + 1
	    	    	    	      #       cat("dfp search intermediate output:","\n")
	    	    	    	      #       cat("f0: ",f0,"  fmin: ",fmin,"   nfcn: ",nfcn,"\n")
	    	    	    	      #       cat("xt: ",xt,"\n")
	    	    	    	      if (xt[i] > xt.inf) {
                                        cat("parameter ", x$label[i], " near upper boundary","\n")}
	    	    	    	      if (xt[i] < -xt.inf) { 
	    	    	    	        cat("parameter ", x$label[i], " near lower boundary","\n")}
	    	    	    	      if (f0 == "NaN") {
	    	    	    	       warning("f0 is NaN")
	    	    	    	       nf <- 10
	    	    	    	       break
	    	    	    	      }
	    	    	    	      if (f0 <= fmin & abs(xt[i]) < xt.inf) {
	    	    	    	       # success	
	    	    	    	       xtf <- xt[i]
	    	    	    	       d <- 3 * d
	    	    	    	       fmin <- f0
	    	    	    	       kg <- 1
	    	    	    	       ns <- ns + 1
	    	    	    	      }
	    	    	    	      else {
	    	    	    	       if (kg == 1) {
	    	    	    	        ns <- 0
	    	    	    	        nf <- 0
	    	    	    	        # failure	
	    	    	    	        break
	    	    	    	       }
	    	    	    	       else {
	    	    	    	        kg <- -1
	    	    	    	        nf <- nf + 1
	    	    	    	        d <- (-0.4) * d
	    	    	    	       }
	    	    	    	      }
	    	    	    	     }
	    	    	    	     if (nf == 10) {
	    	    	    	      d <- 1000 * d
	    	    	    	      xtf <- xbeg
	    	    	    	      g2[i] <- 1
	    	    	    	      negg2 <- negg2 - 1
	    	    	    	     }
	    	    	    	     if (ns == 10) {
	    	    	    	      if (fmin >= fs) {
	    	    	    	       d <- 0.001 * d
	    	    	    	       xtf <- xbeg
	    	    	    	       g2[i] <- 1
	    	    	    	       negg2 <- negg2 - 1
	    	    	    	      }
	    	    	    	     }
	    	    	    	     xt[i] <- xtf
	    	    	    	     dirin[i] <- 0.1 * d
	    	    	    	     fs <- fmin
	    	    	    	    }
	    	    	    }
	    	    }
	    	    ####  provide output
	    	    if (ntry > 6) {
	    	    	    warning("could not find pos. def. covariance - procede with DFP")}
	    	    ntry <- 0
	    	    matgd <- 1
	    	    ####  diagonal matrix
	    	    ####  get sigma and set up loop
	    	    if (isw2 <= 1) {
	    	    	    ntry <- 1
	    	    	    matgd <- 0
	    	    	    v <- matrix(0, npar, npar)
	    	    	    diag(v) <- 2/g2
	    	    }
	    	    if (!all(diag(v) > 0)) {
	    	    	    ntry <- 1
	    	    	    matgd <- 0
	    	    	    v <- matrix(0, npar, npar)
	    	    	    # check whether always g2 > 0 ?
	    	    	    diag(v) <- 2/g2
	    	    }
	    	    xxs <- xt
	    	    sigma <- as.vector(gs %*% (v %*% gs)) * 0.5
	    	    if (sigma < 0) {
	    	    	    cat("covariance matrix is not positive-definite", 
	    	    	    	    "\n")
	    	    	    if (ntry == 0) {
	    	    	    	    #try one more time (setting ntry=1)
	    	    	    	    ntry <- 1
	    	    	    	    matgd <- 0
	    	    	    	    v <- matrix(0, npar, npar)
	    	    	    	    diag(v) <- 2/g2
	    	    	    	    xxs <- xt
	    	    	    	    sigma <- as.vector(gs %*% (v %*% gs)) * 0.5
	    	    	    	    ####  provide output
	    	    	    	    if (sigma < 0) {
	    	    	    	     isw2 <- 0
	    	    	    	     warning("could not find pos. def. covariance")
	    	    	    	     x$est <- btrf(xt, x$low, x$upp)
	    	    	    	     return(list(fmin = fmin, x = x))
	    	    	    	    }
	    	    	    }
	    	    }
	    	    else {
	    	    	    isw2 <- 1
	    	    	    iter <- 0
	    	    }
	    	    x$est <- btrf(xt, x$low, x$upp)
	    	    cat("dfp search intermediate output:", "\n")
	    	    cat("iter: ", iter, "  fmin: ", fmin, "   nfcn: ", 
	    	    	    nfcn, "\n")
	    	    ####  start main loop #######################################
	    	    if(dfp.loop == 0) {cat("start main loop:", "\n")} else {
                      cat("restart main loop:", "\n")}
	    	    f.main <- 0
	    	    #exit main if: f.main > 0
	    	    while (f.main == 0) {
	    	    	    #vector
	    	    	    dirin <- -0.5 * as.vector(v %*% gs)
	    	    	    #scalar
	    	    	    gdel <- as.vector(dirin %*% gs)
	    	    	    ####  linear search along -vg  . . .
	    	    	    xt <- xxs + dirin
	    	    	    xt[xt > xt.inf] <- xt.inf
	    	    	    xt[xt < -xt.inf] <- -xt.inf
	    	    	    f0 <- f(btrf(xt, x$low, x$upp))
	    	    	    nfcn <- nfcn + 1
	    	    	    ####  cat(format(nfcn),"   f0: ",f0,"   ",format(xt),"\n","\n")
	    	    	    ####  change to output on orig. scale
	    	    	    ####  quadr interp using slope gdel
	    	    	    denom <- 2 * (f0 - fmin - gdel)
	    	    	    if (denom <= 0) {
	    	    	    	    slam <- slamax
	    	    	    }
	    	    	    else {
	    	    	    	    slam <- -gdel/denom
	    	    	    	    if (slam > slamax) {
	    	    	    	     ####  cat("slam: ",format(slam),"\n")
	    	    	    	     slam <- slamax
	    	    	    	    }
	    	    	    	    else {
	    	    	    	     if (slam < slamin) 
	    	    	    	      slam <- slamin
	    	    	    	    }
	    	    	    }
	    	    	    # if (abs(slam-1.)  >=  0.1)
	    	    	    if (abs(slam - 1) >= 0.1) {
	    	    	    	    # else go to 70
	    	    	    	    xt <- xxs + slam * dirin
	    	    	    	    xt[xt > xt.inf] <- xt.inf
	    	    	    	    xt[xt < -xt.inf] <- -xt.inf
	    	    	    	    f2 <- f(btrf(xt, x$low, x$upp))
	    	    	    	    nfcn <- nfcn + 1
	    	    	    	    ####   cat(format(nfcn),"  f2: ",f2,"   ",format(xt),"\n")
	    	    	    	    ####   quadr interp using 3 points
	    	    	    	    aa <- fs/slam
	    	    	    	    bb <- f0/(1 - slam)
	    	    	    	    cc <- f2/(slam * (slam - 1))
	    	    	    	    denom <- 2 * (aa + bb + cc)
	    	    	    	    if (denom <= 0) {
	    	    	    	     tlam <- tlamax
	    	    	    	    }
	    	    	    	    else {
	    	    	    	     tlam <- (aa * (slam + 1) + bb * slam + cc)/denom
	    	    	    	     ####  cat("tlam: ",format(tlam),"\n"); cat(f0,f2,fs,'\n')
	    	    	    	     if (tlam > tlamax) {
	    	    	    	      tlam <- tlamax
	    	    	    	     }
	    	    	    	     else {
	    	    	    	      if (tlam < tlamin) 
	    	    	    	       tlam <- tlamin
	    	    	    	     }
	    	    	    	    }
	    	    	    	    xt <- xxs + tlam * dirin
	    	    	    	    xt[xt > xt.inf] <- xt.inf
	    	    	    	    xt[xt < -xt.inf] <- -xt.inf
	    	    	    	    f3 <- f(btrf(xt, x$low, x$upp))
	    	    	    	    nfcn <- nfcn + 1
	    	    	    	    if (f0 >= fmin & f2 >= fmin & f3 >= fmin) {
	    	    	    	     f.main <- 200
	    	    	    	     cat("break 200", "\n")
	    	    	    	     ####	cat(format(nfcn),"  f3: ",f3,"   ",format(xt),"\n")
	    	    	    	     break
	    	    	    	    }
	    	    	    	    if (f0 < f2 & f0 < f3) {
	    	    	    	     slam <- 1
	    	    	    	    }
	    	    	    	    else {
	    	    	    	     if (f2 < f3) {
	    	    	    	      f0 <- f2
	    	    	    	     }
	    	    	    	     else {
	    	    	    	      ## 55?	
	    	    	    	      f0 <- f3
	    	    	    	      slam <- tlam
	    	    	    	     }
	    	    	    	    }
	    	    	    	    dirin <- dirin * slam
	    	    	    	    xt <- xxs + dirin
	    	    	    	    xt[xt > xt.inf] <- xt.inf
	    	    	    	    xt[xt < -xt.inf] <- -xt.inf
	    	    	    }
	    	    	    fmin <- f0
	    	    	    isw2 <- 2
	    	    	    if (sigma + fs - fmin < rostop) {
	    	    	    	    f.main <- 170
	    	    	    	    ####  stop/convergence criteria
	    	    	    	    break
	    	    	    }
	    	    	    apsi.c <- sigma + rho2 + fs - fmin
	    	    	    if (apsi.c <= apsi) {
	    	    	    	    if (trace < vtest) {
	    	    	    	     f.main <- 170
	    	    	    	     break
	    	    	    	    }
	    	    	    }
	    	    	    #non-convergence
	    	    	    if (nfcn - npfn >= nfcnmx) {
	    	    	    	    f.main <- 190
	    	    	    	    break
	    	    	    }
	    	    	    iter <- iter + 1
	    	    	    ####  get gradient and sigma
	    	    	    ####  compute first and second (diagonal) derivatives 
	    	    	    fmin <- f(btrf(xt, x$low, x$upp))
	    	    	    nfcn <- nfcn + 1
	    	    	    ####      cat(format(nfcn),"  ",fmin,"   ",format(xt),"\n")
	    	    	    ####      cat("___________________________________________","\n")
	    	    	    for (i in 1:npar) {
	    	    	    	    vii <- v[i, i]
	    	    	    	    if (vii > 0) {
	    	    	    	     d <- 0.002 * sqrt(abs(vii) * up)
	    	    	    	    }
	    	    	    	    else {
	    	    	    	     d <- 0.02
	    	    	    	    }
	    	    	    	    xtf <- xt[i]
	    	    	    	    xt[i] <- xtf + d
	    	    	    	    fs1 <- f(btrf(xt, x$low, x$upp))
	    	    	    	    nfcn <- nfcn + 1
	    	    	    	    xt[i] <- xtf - d
	    	    	    	    fs2 <- f(btrf(xt, x$low, x$upp))
	    	    	    	    nfcn <- nfcn + 1
	    	    	    	    xt[i] <- xtf
	    	    	    	    g[i] <- (fs1 - fs2)/(2 * d)
	    	    	    	    g2[i] <- (fs1 + fs2 - 2 * fmin)/d^2
	    	    	    }
	    	    	    rho2 <- sigma
	    	    	    vg <- 0.5 * as.vector(v %*% (g - gs))
	    	    	    gvg <- as.vector((g - gs) %*% vg)
	    	    	    delgam <- as.vector(dirin %*% (g - gs))
	    	    	    sigma <- 0.5 * as.vector(g %*% (v %*% g))
	    	    	    if (sigma >= 0) {
	    	    	    	    if (gvg <= 0 | delgam <= 0) {
	    	    	    	     {
	    	    	    	      if (sigma < 0.1 * rostop) {
	    	    	    	       f.main <- 170
	    	    	    	       break
	    	    	    	      }
	    	    	    	      else {
	    	    	    	       f.main <- 1
	    	    	    	       break
	    	    	    	      }
	    	    	    	     }
	    	    	    	    }
	    	    	    }
	    	    	    else {
	    	    	    	    f.main <- 1
	    	    	    	    break
	    	    	    }
	    	    	    ####  update covariance matrix
	    	    	    trace <- 0
	    	    	    vii <- diag(v)
	    	    	    for (i in 1:npar) {
	    	    	    	    for (j in 1:npar) {
	    	    	    	     d <- dirin[i] * dirin[j]/delgam - vg[i] * 
	    	    	    	      vg[j]/gvg
	    	    	    	     v[i, j] <- v[i, j] + 2 * d
	    	    	    	    }
	    	    	    }
	    	    	    if (delgam > gvg) {
	    	    	    	    flnu <- dirin/delgam - vg/gvg
	    	    	    	    for (i in 1:npar) {
	    	    	    	     for (j in 1:npar) {
	    	    	    	      v[i, j] <- v[i, j] + 2 * gvg * flnu[i] * 
	    	    	    	       flnu[j]
	    	    	    	     }
	    	    	    	    }
	    	    	    }
	    	    	    xxs <- xt
	    	    	    gs <- g
	    	    	    trace <- sum(((diag(v) - vii)/(diag(v) + vii))^2)
	    	    	    trace <- sqrt(trace/npar)
	    	    	    fs <- f0
                          } # end main loop ###########################################
	    	    if (f.main == 1) {
	    	    	    x$est <- btrf(xt, x$low, x$upp)
	    	    	    dfp.loop <- dfp.loop + 1
	    	    	    cat("dfp loop: ", dfp.loop, "\n")
	    	    }
	    	    if (f.main == 190) {
	    	    	    #exceeds nfcn
	    	    	    isw1 <- 1
	    	    	    f.main <- 230
                            status <- 1
	    	    	    break
	    	    }
	    	    if (f.main == 200) {
	    	    	    cat("dfp search fails to converge, restart ...", 
	    	    	    	    "\n")
	    	    	    xt <- xxs
	    	    	    x$est <- btrf(xt, x$low, x$upp)
	    	    	    isw2 <- 1
	    	    	    cat("sigma: ", sigma, "\n")
	    	    	    if (sigma < rostop) {
	    	    	    	    f.main <- 170
	    	    	    }
	    	    	    else {
	    	    	    	    if (matgd >= 0) {
	    	    	    	     ###################################
	    	    	    	     dfp.loop <- dfp.loop + 1
	    	    	    	     cat("dfp loop: ", dfp.loop, "\n")
	    	    	    	    }
	    	    	    	    else {
                                      status <- 1
	    	    	    	      break
	    	    	    	    }
	    	    	    }
                     }
	    	    ####  CONVERGENCE
	    	    if (f.main == 170) {
                      status <- 0
	    	    	    cat("DFP search completed with status ",trunc(status),"\n","\n")
	    	    	    isw2 <- 3
	    	    	    if (matgd == 0) {
	    	    	    	    npargd <- npar * (npar + 5)/2
	    	    	    	    if (nfcn - npfn < npargd) {
	    	    	    cat("covariance matrix inaccurate - calculate Hessian","\n")
	    	    	    	     if (isw2 >= 2) {
	    	    	    	      isw2 <- 3
	    	    	    	      cat("perhaps dfp started near minimum - try newton-raphson", "\n")
                                    }
                          }
                                  }
	    	    	    break
	    	    }
                  }
	    fmin <- f(btrf(xt, x$low, x$upp)); nfcn <- nfcn + 1

	    x$est <- btrf(xt, x$low, x$upp)
            ####  compute error (logit scale)
            del <- dqstep(x,f,sens=.01)
            h <- logit.hessian(x,f,del,dapprox=FALSE,nfcn); nfcn <- h$nfcn
            v <- solve(h$ddf)
            
            xtl <- xt-1.96*sqrt(diag(v))
            xtu <- xt+1.96*sqrt(diag(v))
            
            ####  dfp search final output:
            xl <- btrf(xtl, x$low, x$upp)
            xu <- btrf(xtu, x$low, x$upp)
            
	    cat("Optimization Result:", "\n")
	    cat("iter: ", iter, "  fmin: ", fmin, "   nfcn: ", nfcn, 
	    	    "\n")
	    cat("\n")
	    m.out <- cbind(x$label, signif(x$est,8), signif(xl,8), signif(xu,8))
	    dimnames(m.out) <- list(1:npar, c("label", "estimate", "low", "high"))
	    print(m.out, quote = FALSE)
	    cat("\n")
	    return(list(fmin = fmin, label = x$label, est = x$est, status=status, nfcn=nfcn))
}




"dqstep" <-
function(x,f,sens) {
  # fix nfcn counter later
  npar  <- length(x$est)
  step <- .001; dsteps <- rep(step,npar)

  # OBTAIN REFERENCE POINT VALUE
  xt <- ftrf(x$est, x$low, x$upp)
  f0 <- f(x$est)
  
  for(i in 1:npar) {
    stepi <- step
    flag <- 0
    xt.new <- xt
    while(flag==0) {
      flag <- 1
      xt.new[i] <- xt[i]-stepi; x1 <- -stepi
      f1 <- f(btrf(xt.new, x$low, x$upp))
      xt.new[i] <- xt[i]+stepi; x2 <- stepi
      f2 <- f(btrf(xt.new, x$low, x$upp))

      # handle exceptions
      if (is.na(f2)) f2 <- Inf
      if(f2==Inf | f2==-Inf) {
        warning('Infs - reducing step size')
        stepi <- stepi/10; flag <- 0}
      if(f2==f0 & f1==f0) {
        cat('increasing step size','\n')
        stepi <- stepi*10; flag <- 0}
      if(abs(f2-f0) > (.5 * sens) | abs(f1-f0) > (.5 * sens)) {
        cat('reducing step size','\n')
        stepi <- stepi/10; flag <- 0
        # cat(stepi,f1,f2,'\n')
      }
    }
    
    b <- ((f1-f0)*x2-(f2-f0)*x1)/(x1*x1*x2-x2*x2*x1)
    a <- (f1-f0)/x1-b*x1

    # ***roots
    r <- a*a+4*b*sens
    if(r < 0 | is.na(r) | b ==0) {
      warning('oops: unable to find stepsize, use default')
      cat('problem with ',x$label[i],'\n')
      break
    }
      
    xs1 <- 0.5*(-a-sqrt(a*a+4*b*sens))/b
    xs2 <- 0.5*(-a+sqrt(a*a+4*b*sens))/b

    if(abs(xs1) <= abs(xs2)) {xs <- xs1} else { xs <- xs2}

    # *** see where we end up
    xt.new[i] <- xt[i]+xs
    f2 <- f(btrf(xt.new, x$low, x$upp))
    if (is.na(f2)) f2 <- Inf
    if(f2==Inf | f2==-Inf) {
      warning('oops: unable to find stepsize, use default')
      break
    }
    dsteps[i] <- xs
    # cat('DSTEP:',x$label[i],signif(xs,6),(f2-f0),'\n')
  }
  return(dsteps)
}

"ftrf" <-
function(x,xl,xu) {

####  forward transformation
####  this assumes logit transformations of the parameters
####  bounded from below by xl and from above by xu
  if(any((x-xl) <= 0)) stop('ftrf requires x > xl')
  if(any((xu-x) <= 0)) stop('ftrf requires x < xu')
  return(log((x-xl)/(xu-x)))
}
"logit.hessian" <-
function(x=x,f=f,del=rep(.002,length(x$est)),dapprox=FALSE,nfcn=0) {
  # computes the hessian of f on logit scale
  small <- 1.e-8
  npar <- length(x$est)
  nlm <- 2*npar*npar
  np2 <- 2*npar
  
  xt <- ftrf(x$est, x$low, x$upp)
  f0 <- f(x$est); nfcn <- nfcn + 1
  # cat('compute internal Hessian with function value at:',f0,'\n')

  # *** INITIALIZE NEWTON BY CALCULATION OF A
  #     NLM POINT SIMPLEX GEOMETRY IN N DIM PARAMETER SPACE
  xn <- matrix(rep(xt,nlm),ncol=nlm)
  
  # ***   ON AXES
  # ***	compute delta distance in each direction

  for(i in 1:npar) {
  j.even <- 2*i
  j.odd  <- j.even-1
  xn[i,j.even] <- xn[i,j.even]-del[i]
  xn[i,j.odd]  <- xn[i,j.odd] +del[i]
  }

  if(dapprox==FALSE) {
   # ***   OFF AXIS
   if(npar > 1) {
     mc <- np2+1
     for(i in 2:npar) {
       for(j in 1:(i-1)) {
         xn[i,mc] <- xn[i,mc]+del[i]
         xn[j,mc] <- xn[j,mc]+del[j]
         mc <- mc+1
         xn[i,mc] <- xn[i,mc]-del[i]
         xn[j,mc] <- xn[j,mc]+del[j]
         mc <- mc+1
         xn[i,mc] <- xn[i,mc]-del[i]
         xn[j,mc] <- xn[j,mc]-del[j]
         mc <- mc+1
         xn[i,mc] <- xn[i,mc]+del[i]
         xn[j,mc] <- xn[j,mc]-del[j]
         mc <- mc+1
       }}
   }
   
   f.vec <- numeric(nlm) 
   for(i in 1:nlm) {
     f.vec[i] <- f(btrf(xn[,i], x$low, x$upp))
     nfcn <- nfcn + 1
   }
 } else { # ignore off diagonal elements

   f.vec <- numeric(np2) 
   for(i in 1:np2) {
     f.vec[i] <- f(btrf(xn[,i], x$low, x$upp))
     nfcn <- nfcn + 1
   }
 }
  
  # ***   FIRST AND DIAGONAL SECOND DERIVATIVES
 
  i <- 1:npar; i.even <- 2*i; i.odd <- i.even-1
 
  df <- (f.vec[i.even]-f.vec[i.odd])/2/del
  if(npar > 1) {ddf <- diag((f.vec[i.even]+f.vec[i.odd]-2*f0)/(del**2))/2 }
  else {ddf <- (f.vec[i.even]+f.vec[i.odd]-2*f0)/(del**2)/2 }  
  # print(format(ddf),quote=FALSE)
 
  if(dapprox==FALSE) {
  # ***   SECOND DERIVATIVES
 
  if(npar > 1) {
 
    mc <- np2+1
    for(i in 2:npar) {
      for(j in 1:(i-1)) {
        mc1 <- mc; mc2 <- mc1+1; mc3 <- mc2+1; mc4 <- mc3+1
        ddf[i,j] <- (f.vec[mc1]+f.vec[mc3]-f.vec[mc2]-f.vec[mc4])/del[i]/del[j]/4
        mc <- mc4+1
      }}
    ddf <- ddf+t(ddf)
  }
} else {ddf <- 2 * ddf}
 
  # ***   EIGEN VALUES and MATRIX INVERSION
 
  eig <- eigen(ddf)
  # cat('eigen values: ',format(sort(eig$values)),'\n')
  if(any(eig$values) < 0) {
    warning('hessian not pos. definite')
  }
  if(any(abs(eig$values)) < small) {
    warning('hessian may be singular')
  }
  # *** ADJUSTMENT OF del (feature not implemented) 
  # ddf.diag <- diag(ddf)
  # del[ddf.diag > 0] <- .002/sqrt(ddf.diag[ddf.diag > 0])

  return(list(df=df,ddf=ddf,nfcn=nfcn,eigen=eig$values))
}
"mcmc" <-
function (x, nlogf, m1, m2=min(m1,m3), m3, scl1=0.5, scl2=2, skip=1, covm=0, nfcn = 0, plot=FALSE)
{
	    #     MCMC/MH sampler for R
	    #     This module is part of the Bhat likelihood exploration tool.

	    #     This program is free software; you can redistribute it and/or modify
	    #     it under the terms of the GNU General Public License as published by
	    #     the Free Software Foundation; either version 2 of the License, or
	    #     (at your option) any later version.
	    #     This program is distributed in the hope that it will be useful,
	    #     but WITHOUT ANY WARRANTY; without even the implied warranty of
	    #     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
	    #     GNU General Public License for more details.

            #     E. Georg Luebeck (gluebeck@fhcrc.org)

	    if (!is.list(x)) {
	    	    cat("x is not a list! see help file", "\n")
	    	    return()
	    }
	    names(x)[1] <- "label"
	    names(x)[2] <- "est"
	    names(x)[3] <- "low"
	    names(x)[4] <- "upp"
	    npar <- length(x$est)
            
	    ####  objects:
	    if (npar <= 0) stop('no. of parameters must be >= 1')

            small <- 1.e-8
            sens <- 0.1
            xinf <- 20.

            # initialize graphical output
            if(plot==TRUE) {
            par(mfrow=c(npar+1,1),bg="grey")}

            # *** initialize f.old
	    cat(date(), "\n")
	    xt <- ftrf(x$est, x$low, x$upp)
	    f.old <- nlogf(x$est); nfcn <- nfcn + 1
            f.first <- f.old

            acc.count <- 0

            # *** if COVM is not provided
            if(!is.matrix(covm)) {

            # initialize monitors
            f.mon <- numeric(m1+m3); x.mon <- matrix(0,ncol=npar,nrow=m1+m3)

            cat('trying a proposal COVM using the Hessian','\n')
            del <- dqstep(x,nlogf,sens)
            h <- logit.hessian(x,nlogf,del,dapprox=FALSE,nfcn)  # returns list: $df,$ddf,$nfcn
            nfcn <- h$nfcn
            ddf <- h$ddf

            # ***   EIGEN VALUES and MATRIX INVERSION
            eig <- eigen(ddf)
            
            if(any(eig$values < small)) {
                  cat('Hessian may not be positive definite','\n')
                  cat('trying a proposal COVM using dqstep','\n')
                  del <- dqstep(x,nlogf,4.)
                  eig <- .1/del/del; ddf <- diag(eig); eig <- eigen(ddf)
                } 

            # ggf <- 0.5 * solve(ddf,diag(1,npar),tol=1.e-10)
            
            # cat('unity test:','\n'); print(format(ggf %*% ddf),quote=FALSE)
            cat('first pilot chain: using MLEs and log-transformed covariance','\n','\n')
            
            nc <- 0
            for (n in 1:m1) {

              accept <- 1

              # *** compute proposal x' (=yt). If unacceptable, reject
              dx <-  eig$vectors %*% rnorm(npar,0,1/(0.5*sqrt(eig$values)))
              yt <- xt + scl1 * dx
                                          
              # *** get log-likelihood
              f.new <- nlogf(btrf(yt, x$low, x$upp)); nfcn <- nfcn + 1
                                      
              # *** boundary checks from within func ...

              # *** R ratio 
              if(any(abs(yt) > xinf)) {
                cat('cycle',n,': mcmc close to boundary, move rejected!','\n'); accept <- 0} else {
              accept <- min(accept,exp(-f.new+f.old))}
              
              if(accept == 1) {xt <- yt; f.old <- f.new; acc.count <- acc.count+1
                             xt.maxl <- yt; f.maxl <- f.new # approximate/search max-likelihood
                             } else {
                               if(runif(1) <= accept) {
                                         xt <- yt; f.old <- f.new; acc.count <- acc.count+1}
                             }

              f.mon[n] <- f.old; x.mon[n,] <- btrf(xt, x$low, x$upp)

              # *** regular output and graphical monitor
              if(n%%(100*skip) == 0) {
                m.out <- c("n:",signif(n,3),"acceptance:",signif(acc.count/100/skip,3),"-log lkh:",signif(f.new,8),signif(x.mon[n,],6))
                cat(m.out,'\n')
                acc.count <- 0
                n.skip <- seq(skip,n,skip)

                if(plot==TRUE) {
                par(mar=c(0, 5, 0.1, 4) + 0.1)
                plot(f.mon[n.skip], type='l', xlab = " ", ylab = "-log-density",col=2)
                for (i in 1:(npar-1)) {
                  par(mar=c(0, 5, 0, 4) + 0.1)
                  plot(x.mon[n.skip,i], type='l', xlab = " ", ylab = x$label[i], col=2)
                }
                par(mar=c(0.1, 5, 0, 4) + 0.1)
                plot(x.mon[n.skip,npar], type='l', xlab = " ", xaxt='n', ylab = x$label[npar], col=2)
                nc <- nc + 1
              }
              }
            }

            # *** obtain empirical covariance of increments
            covm <- cov(x.mon[2:m1,]-x.mon[1:(m1-1),]) # now redundant
            cat('\n','\n')
            # print(covm,quote=FALSE)
            # cat('\n','\n')

            # *** pilot run 2 *** includes m2 cycles for incremental adjustment of covm
            cat('second pilot run:','\n','\n')
            eig <- eigen(covm)
            cat('eigen values:',eig$values,'\n','\n')

          } else {
            # covm for mvn proposal distribution given as input  
            if(nrow(covm)!=length(x$est)) {stop('dimension of covm not specified correctly')}
            nc <- 0 
            # re-initialize monitors
            m1 <- 1
            f.mon <- numeric(m1+m3); x.mon <- matrix(0,ncol=npar,nrow=m1+m3)
            eig <- eigen(covm); x.mon[1,] <- x$est; f.mon <- f.first
          }
            acc.count <- 0
            for (n in (m1+1):(m1+m3)) {

              x.mon[n,] <- x.mon[n-1,]; f.mon[n] <- f.mon[n-1]
              accept <- 1

              # *** compute proposal x' (=yt). If unacceptable, reject
              dx <-  eig$vectors %*% rnorm(npar,0,sqrt(eig$values))
              y <- x.mon[n-1,] + scl2 * dx
              if(any((y-x$low) < small)) {cat('mcmc move rejected (lower bound)','\n'); accept <- 0}
              if(any((x$upp-y) < small)) {cat('mcmc move rejected (upper bound)','\n'); accept <- 0}

              # *** boundary checks from within func ...

              if(accept > 0) {
              # *** get log-likelihood
              f.new <- nlogf(y); nfcn <- nfcn + 1
                                      
              # *** acceptance ratio and updates
              accept <- min(accept,exp(-f.new+f.mon[n]))
              if(accept == 1) {x.mon[n,] <- y; f.mon[n] <- f.new; acc.count <- acc.count+1
                             x.maxl <- y; f.maxl <- f.new # approximate/search max-likelihood
                             } else {
                               if(runif(1) <= accept) {
                                 x.mon[n,] <- y; f.mon[n] <- f.new; acc.count <- acc.count+1
                               }
                             }
            }

              # *** regular output and graphical monitor
              if(n%%(100*skip) == 0) {
                m.out <- c("n:",signif(n,3),"acceptance:",signif(acc.count/100/skip,3),"-log lkh:",signif(f.new,8),signif(x.mon[n,],6))
                cat(m.out,'\n')
                acc.count <- 0

                n.skip.1 <- seq(skip,n,skip)
                n.skip.2 <- seq(skip,min(n,m1+m2),skip)

                if(plot==TRUE) {
                if(m1 > 1) {brncol <- 3} else {brncol <- 2}
                  
                par(mar=c(0, 5, 0.1, 4) + 0.1)
                plot(f.mon[n.skip.1], type='l', xlab = " ", ylab = "-log-density",col=2)
                lines(f.mon[n.skip.2],col=brncol) #pilot cycles
                for (i in 1:(npar-1)) {
                  par(mar=c(0, 5, 0, 4) + 0.1)
                  plot(x.mon[n.skip.1,i], type='l', xlab = " ", ylab = x$label[i], col=2)
                  lines(x.mon[n.skip.2,i],col=brncol) #pilot cycles
                }
                par(mar=c(0.1, 5, 0, 4) + 0.1)
                plot(x.mon[n.skip.1,npar], type='l', xlab = " ", xaxt='n', ylab = x$label[npar], col=2)
                lines(x.mon[n.skip.2,npar],col=brncol) #pilot cycles
                
                nc <- nc + 1
              }

                if(m1 > 1) { #note: when covm is passed to mcmc, m1 is set to 1
                # update covariance using sampled increments (m1+1):n
                if(n <= m1+m2) {
                covm <- cov(x.mon[2:n,]-x.mon[1:(n-1),])
                eig <- eigen(covm)
                if(any(eig$values < small)) {warning('covm nearly singular')}
              }
              }
              }
            }
            return(list(f=f.mon,mcmc=x.mon,covm=covm))
            }











"newton" <-
function (x, f, eps=1e-1, itmax=10, relax=0, nfcn = 0) 
{
	    #     General Newton-Raphson module written in R. 
	    #     This module is part of the Bhat likelihood exploration tool.

	    #     This program is free software; you can redistribute it and/or modify
	    #     it under the terms of the GNU General Public License as published by
	    #     the Free Software Foundation; either version 2 of the License, or
	    #     (at your option) any later version.
	    #     This program is distributed in the hope that it will be useful,
	    #     but WITHOUT ANY WARRANTY; without even the implied warranty of
	    #     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
	    #     GNU General Public License for more details.

            #     E. Georg Luebeck, Fred Hutchinson Cancer Research Center
	    #     Fairview Avenue N, Seattle, WA 98109-1024

	    if (!is.list(x)) {
	    	    cat("x is not a list! see help file", "\n")
	    	    return()
	    }
	    names(x)[1] <- "label"
	    names(x)[2] <- "est"
	    names(x)[3] <- "low"
	    names(x)[4] <- "upp"
	    npar <- length(x$est
                           )
	    ####  objects:
	    if (npar <= 0) {
	    	    warning("no. of parameters < 1")
	    	    exit
	    }

            small <- 1.e-8
            nlm <- (npar*npar+3*npar)/2
            np2 <- 2*npar
            sens <- 0.1
            flag.hessian <- 1

	    ####  first call 
	    cat(date(), "\n")
	    xt <- ftrf(x$est, x$low, x$upp)
	    f0 <- f(x$est); nfcn <- nfcn + 1

            for(n in 1:itmax) {

              cat(nfcn, '  fmin: ', format(f0), "   ", format(x$est), '\n')

                # *** INITIALIZE NEWTON BY CALCULATION OF A
                #     NLM POINT SIMPLEX GEOMETRY IN NPAR DIM PARAMETER SPACE
                xn <- matrix(rep(xt,nlm),ncol=nlm)
                # ***   ON AXES
                if(n==1) del <- dqstep(x,f,sens)
                if(flag.hessian==1) {
                  h <- logit.hessian(x,f,del,dapprox=FALSE,nfcn)  # returns list: $df,$ddf,$nfcn
                  nfcn <- h$nfcn
                  ddf <- h$ddf; df <- h$df
                } else {

                # ***	compute delta distance in each direction
                for(i in 1:npar) {
                  j.even <- 2*i
                  j.odd  <- j.even-1
                  xn[i,j.even] <- xn[i,j.even]-del[i]
                  xn[i,j.odd]  <- xn[i,j.odd] +del[i]
                }

                # ***   OFF AXIS
                if(npar > 1) {
                  mc <- np2+1
                  for(i in 2:npar) {
                    for(j in 1:(i-1)) {
                      xn[i,mc] <- xn[i,mc]+del[i]
                      xn[j,mc] <- xn[j,mc]+del[j]
                      mc <- mc+1
                    }}
                }

                f.vec <- numeric(nlm) 
                for(i in 1:nlm) {
                  f.vec[i] <- f(btrf(xn[,i], x$low, x$upp))
	    	  nfcn <- nfcn + 1
                }
                
                # ***   FIRST AND DIAGONAL SECOND DERIVATIVES
                i <- 1:npar
                i.even <- 2*i
                i.odd <- i.even-1

                df <- (f.vec[i.even]-f.vec[i.odd])/2/del
                ddf <- diag((f.vec[i.even]+f.vec[i.odd]-2*f0)/del/del/2)

                # ***   SECOND DERIVATIVES

                if(npar > 1) {
                  mc <- np2+1
                  for(i in 2:npar) {
                    for(j in 1:(i-1)) {
                      ip <- 2*i-1; jp <- 2*j-1
                      ddf[i,j] <- (f0+f.vec[mc]-f.vec[jp]-f.vec[ip])/del[i]/del[j]
                      mc <- mc+1
                    }}

                  ddf <- ddf+t(ddf)
                }
              
                # ***   EIGEN VALUES and MATRIX INVERSION

                eig <- eigen(ddf)
                cat('approx. eigen values: ',format(eig$values),'\n')
                if(any(eig$values < 0)) {
                  warning('approx. hessian not pos. definite')
                }
                if(any(abs(eig$values) < small)) {
                  warning('approx. hessian may be singular')
                }
              }
                # cat('compare:','\n')
                # print(format(h$ddf),quote=FALSE)
                # print(format(ddf),quote=FALSE)
                      
                b <- diag(1,npar)
                # *** if inversion fails need to return or reinitialize NR (not implemented)
                ggf <- solve(ddf,b,tol=1.e-10)

                # cat('unity test:','\n'); print(format(ggf %*% ddf),quote=FALSE)
              
                # ---	variable Newton stepping

                disc <- ggf  %*% df
                disc0 <- disc
                rneps <- 1; nepsc <- 0; xt0 <- xt; f1 <- f0

                xt <- xt0 + disc
                # cat(xt0,'\n',xt,'\n',df,'\n')
                tmp <- numeric(npar)
                tmp[abs(df) <= eps] <- 1
                nz <- sum(tmp)
                f0 <- f(btrf(xt, x$low, x$upp)); nfcn <- nfcn + 1

                if(relax > 0) {
                while(nepsc < 8 && f0 > f1 ) {
                  rneps <- .75*rneps
                  cat('reducing step size: ',format(rneps,f0),'\n')
                  nepsc <- nepsc+1
                  disc <- rneps * disc0
                  xt <- xt0 + disc
                  f0 <- f(btrf(xt, x$low, x$upp)); nfcn <- nfcn + 1
                 }
                }

                if(nepsc == 8) warning('problem finding lower function value, will continue')
                
                # ---   STOP-CRITERION FOR ITERATION --------------------------------------
                if(nz == npar && min(h$eigen) > small) {

                status <- 'converged'

                # ---	compute Hessian symmetrically:
                x$est <- btrf(xt,x$low,x$upp)  
                # h <- logit.hessian(x,f,del,dapprox=FALSE,nfcn)  # returns list: $df,$ddf,$nfcn
                b <- diag(1,npar)
                ggf <- solve(ddf,b,tol=1.e-10)

                xtl <- rep(NA,npar); xtu <- rep(NA,npar)
                se <- numeric(1:npar)
                ggf.diag <- diag(ggf)
                se[ggf.diag >= 0] <- sqrt(ggf.diag[ggf.diag >= 0])
                xtl <- xt - 1.96 * se
                xtu <- xt + 1.96 * se

                x$est <- btrf(xt, x$low, x$upp)
                xlow  <- btrf(xtl, x$low, x$upp)
                xupp  <- btrf(xtu, x$low, x$upp)

                cat('\n')  
                cat('Bhat run: ',date(),' status: ',status,'\n')
                # cat('Optimization Result:', '\n')
                cat("iter: ", n, "  fmin: ", format(f0), "   nfcn: ", nfcn, "\n")
                cat("\n")
                m.out <- cbind(x$label, round(x$est,6), round(h$df,6), round(xlow,6), round(xupp,6))
                dimnames(m.out) <- list(1:npar, c("label", "estimate", "log deriv", "lower 95%" , "upper 95%"))

                print(m.out, quote=FALSE)
                cat("\n")
                return(list(fmin = f0, label = x$label, est = x$est, low = xlow, upp = xupp))

              }

                # --- graphical diagnostic (do be done)

                # ***	DURING NON-CONVERGEING CYCLES
                status <- 'non-converged'
                xtl <- rep(NA,npar); xtu <- rep(NA,npar)
                se <- numeric(1:npar)
                ggf.diag <- diag(ggf)
                se[ggf.diag >= 0] <- sqrt(ggf.diag[ggf.diag >= 0])
                xtl <- xt - 1.96 * se
                xtu <- xt + 1.96 * se

                x$est <- btrf(xt, x$low, x$upp)
                xlow  <- btrf(xtl, x$low, x$upp)
                xupp  <- btrf(xtu, x$low, x$upp)

                if(n == itmax) {
                  cat('\n')
                  cat('Bhat run: ',date(),' status: ',status,'\n')
                  # cat('Optimization Result:', '\n')
                } 
              
                cat('\n',"iter: ", n, "  fmin: ", f0, "   nfcn: ", nfcn, 
          	    "\n")
                cat("\n")
                vbar <- rep("|",npar)
                m.out <- cbind(x$label, round(x$est,6), round(df,6), round(xlow,6), round(xupp,6))
                dimnames(m.out) <- list(1:npar, c("label", "estimate", "log deriv", "lower 95%" , "upper 95%"))
                print(m.out, quote = FALSE)
                cat("\n")
            }
                return(list(fmin = f0, label = x$label, est = x$est, low = xlow, upp = xupp))
        }
#       Program:     plkhci.R
#                    modified newton-raphson optimizer
#                    to obtain profile-likelihood-based confidence intervals
#
#       Method:      D.J. VENZON and S.H. MOOLGAVKAR
#                    " A method for computing profile-likelihood-based
#                      confidence intervals "
#                    Journal of the Royal Statistical Society, Series C
#                       volume 37, no.1, 1988, pp. 87-94
#                    not yet fully implemented, but may converge
#                    
#       Notes:       p:  identifies which parameter (among npar) is to be used for profile
#                    likelihood calc. 
#                    npar:	total number of parameters in model

plkhci <- function(x, nlogf, label, prob=0.95, eps=.001, nmax=10, nfcn = 0) {
  
  if (!is.list(x)) {
  	    cat("x is not a list! see help file", "\n")
  	    return()
  }
  names(x)[1] <- "label"
  names(x)[2] <- "est"
  names(x)[3] <- "low"
  names(x)[4] <- "upp"
  npar <- length(x$est)

  if (npar <= 0) {
  	    warning("no. of parameters < 1")
  	    exit
  }

  small <- 1.e-8
  nd1 <- npar-1    
  sens <- .5
  rhs <- numeric(npar)
  disc <- numeric(npar)
  qchi <- qchisq(prob,1)

### check selected parameter
  if(!is.character(label)) stop('label needs to be of type character')
  p <- match(label,x$label)
  if(is.na(p)) stop('label not found')
  
  xt <- ftrf(x$est, x$low, x$upp)
  fmle <- nlogf(x$est); nfcn <- nfcn + 1

  cat(' neg. log. likelihood: ',fmle,'\n','\n')
  cat(' will atempt to compute both bounds (+/- direction)','\n')

        del <- dqstep(x,nlogf,sens)

        # *** shoot along tangent?
        # *** initialize first step

        h <- logit.hessian(x,nlogf,del,dapprox=FALSE,nfcn)  # returns list: $df,$ddf,$nfcn
        nfcn <- h$nfcn
        ddf <- h$ddf; df <- h$df

        # ****  construct ddl
	ddl <- ddf[-p,-p]

        # ****  construct dbo
        dbo <- ddf[p,-p]

        # ****  construct dbb
	dbb0 <- ddf[p,p]

        # ***   matrix inversion: ddl inverse
        b <- diag(1,nd1)
        # *** if inversion fails need to return or reinitialize NR (not implemented)
        ggl <- solve(ddl,b,tol=1.e-16)
        dbl0 <- ggl %*% dbo
        tmp0 <- t(dbo) %*% dbl0 
        cat('\n')
  
      for(idir in c(1,-1)) {

        xt0 <- xt
        xt0[p] <- xt0[p]  +idir *       sqrt(qchi/(dbb0-tmp0))/2
        xt0[-p] <- xt0[-p]-idir * as.vector(dbl0) * sqrt(qchi/(dbb0-tmp0))/2

        f0 <- nlogf(btrf(xt0, x$low, x$upp)); nfcn <- nfcn + 1
        
        if(idir ==1)  cat('trying lower bound ------------------------', '\n')
        if(idir ==-1) cat('trying upper bound ------------------------', '\n')
        cat('starting at:   ',2*(f0-fmle),'\n')
        cat('initial guess: ',btrf(xt0, x$low, x$upp),'\n')

# *** search via newton-raphson

      cat('\n')
      cat('begin Newton-Raphson search for profile lkh conf. bounds:','\n')
      cat('eps value for stop criterium:',eps,'\n')
      cat('nmax                        :',nmax,'\n')  

        for (n in 1:nmax) {

        x0 <- x; x0$est <- btrf(xt0, x$low, x$upp)
        f0 <- nlogf(x0$est); nfcn <- nfcn + 1
        h <- logit.hessian(x0,nlogf,del,dapprox=FALSE,nfcn)  # returns list: $df,$ddf,$nfcn
        nfcn <- h$nfcn
        ddf <- h$ddf; df <- h$df

        # ****  construct ddl
	ddl <- ddf[-p,-p]

        # ****  construct dbo
        dbo <- ddf[p,-p]

        # ***  construct matrix fdd
        fdd <- ddf 

        fdd[1,1] <- -df[p]
        rhs[1]  <- (fmle - f0 + 0.5*qchi)
        fdd[1,2:npar] <- -df[-p]
        rhs[2:npar] <- df[-p]
        fdd[-1,-1] <- ddl #lower right block   
        fdd[2:npar,1] <- dbo
      
        b <- diag(1,npar)
        ggl <- solve(fdd,b,tol=1.e-16)

        # *** shoot ahead

        if(n <=2) rneps <- .2
        if(n>2 && n <=4) rneps <- .8
        if(n > 4) rneps <- 1.
        
        nepsc <- 0
        xt00 <- xt0
      
        disc0 <- disc <- as.vector(ggl %*% rhs)

        rhs1 <- 1000.
        #while(nepsc < 8 &&  rhs1 > qchi) {
        disc <- rneps * disc0

        nz <- 0
        xt0[p]  <- xt00[p]+disc[1]
        xt0[-p] <- xt00[-p]+disc[2:npar]
                tmp <- numeric(npar)
                tmp[abs(disc) <= eps] <- 1
                nz <- sum(tmp)

        # *** update f0
        f0 <- nlogf(btrf(xt0, x$low, x$upp)); nfcn <- nfcn + 1
        rhs1 <- abs(fmle-f0+0.5*qchi)

        #if(rhs1 > qchi) {
        #   rneps <- .75 * rneps
        #   cat('reducing step size:',rneps,rhs1,'\n') 
        #   nepsc <- nepsc + 1
        # }

        #if(nepsc ==8) cat('problem closing in, will continue.','\n')

        estimate <- btrf(xt0, x$low, x$upp)
        # ***   stop-criterion for iteration
        if(nz ==  npar && nepsc  ==  0) {
        cat('\n','CONVERGENCE: ',trunc(n),' iterations','\n')
        status <- 'converged'
        cat('\n')
        cat('chisquare value is: ',2*(f0-fmle),'\n')
        cat('confidence bound of ',x$label[p],' is ',estimate[p],'\n')
        cat('log derivatives:    ',df[-p],'\n')
        m.out <- cbind(x$label, signif(estimate,6), signif(df,6), signif(diag(ddf),6))
        dimnames(m.out) <- list(1:npar, c("label", "estimate", "log deriv", "log curv"))
        print(m.out, quote = FALSE)
        cat('\n');break
      }
        # cat('\n')
        # cat('chisquare value is: ',2*(f0-fmle),'\n')
        # cat('confidence bound of ',x$label[p],' is ',estimate[p],'\n')
        # cat('log derivatives:    ',df[-p],'\n')
      }
      if(idir==-1) {xlow <- estimate[p]} else {xupp <- estimate[p]}
      }
        return(c(xlow,xupp))
}








